xref: /petsc/src/mat/impls/h2opus/cuda/math2opus.cu (revision cc4c1da905d89950b196b027190941013bd3e15a)
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>
15*cc4c1da9SBarry Smith   #include <petsc/private/deviceimpl.h> /*I "petscmat.h"   I*/
16f17f7ee2SSatish Balay   #include <petscsf.h>
17f17f7ee2SSatish Balay 
18f17f7ee2SSatish Balay /* math2opusutils */
190dd791a8SStefano Zampini PETSC_INTERN PetscErrorCode MatDenseGetH2OpusStridedSF(Mat, PetscSF, PetscSF *);
20f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode VecSign(Vec, Vec);
21f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode VecSetDelta(Vec, PetscInt);
22f17f7ee2SSatish Balay PETSC_INTERN PetscErrorCode MatApproximateNorm_Private(Mat, NormType, PetscInt, PetscReal *);
23f17f7ee2SSatish Balay 
24f17f7ee2SSatish Balay   #define MatH2OpusGetThrustPointer(v) thrust::raw_pointer_cast((v).data())
25f17f7ee2SSatish Balay 
26f17f7ee2SSatish Balay   /* Use GPU only if H2OPUS is configured for GPU */
27f17f7ee2SSatish Balay   #if defined(PETSC_HAVE_CUDA) && defined(H2OPUS_USE_GPU)
28f17f7ee2SSatish Balay     #define PETSC_H2OPUS_USE_GPU
29f17f7ee2SSatish Balay   #endif
30f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
31f17f7ee2SSatish Balay     #define MatH2OpusUpdateIfNeeded(A, B) MatBindToCPU(A, (PetscBool)((A)->boundtocpu || (B)))
32f17f7ee2SSatish Balay   #else
333ba16761SJacob Faibussowitsch     #define MatH2OpusUpdateIfNeeded(A, B) PETSC_SUCCESS
34f17f7ee2SSatish Balay   #endif
35f17f7ee2SSatish Balay 
36f17f7ee2SSatish Balay // TODO H2OPUS:
37f17f7ee2SSatish Balay // DistributedHMatrix
38f17f7ee2SSatish Balay //   unsymmetric ?
39f17f7ee2SSatish Balay //   transpose for distributed_hgemv?
40f17f7ee2SSatish Balay //   clearData()
41f17f7ee2SSatish Balay // Unify interface for sequential and parallel?
42f17f7ee2SSatish Balay // Reuse geometric construction (almost possible, only the unsymmetric case is explicitly handled)
43f17f7ee2SSatish Balay //
449371c9d4SSatish Balay template <class T>
459371c9d4SSatish Balay class PetscPointCloud : public H2OpusDataSet<T> {
46f17f7ee2SSatish Balay private:
47f17f7ee2SSatish Balay   int            dimension;
48f17f7ee2SSatish Balay   size_t         num_points;
49f17f7ee2SSatish Balay   std::vector<T> pts;
50f17f7ee2SSatish Balay 
51f17f7ee2SSatish Balay public:
52d71ae5a4SJacob Faibussowitsch   PetscPointCloud(int dim, size_t num_pts, const T coords[])
53d71ae5a4SJacob Faibussowitsch   {
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 
71d71ae5a4SJacob Faibussowitsch   PetscPointCloud(const PetscPointCloud<T> &other)
72d71ae5a4SJacob Faibussowitsch   {
73f17f7ee2SSatish Balay     size_t N         = other.dimension * other.num_points;
74f17f7ee2SSatish Balay     this->dimension  = other.dimension;
75f17f7ee2SSatish Balay     this->num_points = other.num_points;
76f17f7ee2SSatish Balay     this->pts.resize(N);
779371c9d4SSatish Balay     for (size_t i = 0; i < N; i++) this->pts[i] = other.pts[i];
78f17f7ee2SSatish Balay   }
79f17f7ee2SSatish Balay 
809371c9d4SSatish Balay   int getDimension() const { return dimension; }
81f17f7ee2SSatish Balay 
829371c9d4SSatish Balay   size_t getDataSetSize() const { return num_points; }
83f17f7ee2SSatish Balay 
84d71ae5a4SJacob Faibussowitsch   T getDataPoint(size_t idx, int dim) const
85d71ae5a4SJacob Faibussowitsch   {
86f17f7ee2SSatish Balay     assert(dim < dimension && idx < num_points);
87f17f7ee2SSatish Balay     return pts[idx * dimension + dim];
88f17f7ee2SSatish Balay   }
89f17f7ee2SSatish Balay 
90d71ae5a4SJacob Faibussowitsch   void Print(std::ostream &out = std::cout)
91d71ae5a4SJacob Faibussowitsch   {
92f17f7ee2SSatish Balay     out << "Dimension: " << dimension << std::endl;
93f17f7ee2SSatish Balay     out << "NumPoints: " << num_points << std::endl;
94f17f7ee2SSatish Balay     for (size_t n = 0; n < num_points; n++) {
959371c9d4SSatish Balay       for (int d = 0; d < dimension; d++) out << pts[n * dimension + d] << " ";
96f17f7ee2SSatish Balay       out << std::endl;
97f17f7ee2SSatish Balay     }
98f17f7ee2SSatish Balay   }
99f17f7ee2SSatish Balay };
100f17f7ee2SSatish Balay 
1019371c9d4SSatish Balay template <class T>
1029371c9d4SSatish Balay class PetscFunctionGenerator {
103f17f7ee2SSatish Balay private:
1048434afd1SBarry Smith   MatH2OpusKernelFn *k;
105f17f7ee2SSatish Balay   int                dim;
106f17f7ee2SSatish Balay   void              *ctx;
107f17f7ee2SSatish Balay 
108f17f7ee2SSatish Balay public:
1098434afd1SBarry Smith   PetscFunctionGenerator(MatH2OpusKernelFn *k, int dim, void *ctx)
110d71ae5a4SJacob Faibussowitsch   {
1119371c9d4SSatish Balay     this->k   = k;
1129371c9d4SSatish Balay     this->dim = dim;
1139371c9d4SSatish Balay     this->ctx = ctx;
114f17f7ee2SSatish Balay   }
115d71ae5a4SJacob Faibussowitsch   PetscFunctionGenerator(PetscFunctionGenerator &other)
116d71ae5a4SJacob Faibussowitsch   {
1179371c9d4SSatish Balay     this->k   = other.k;
1189371c9d4SSatish Balay     this->dim = other.dim;
1199371c9d4SSatish Balay     this->ctx = other.ctx;
1209371c9d4SSatish Balay   }
1219371c9d4SSatish Balay   T operator()(PetscReal *pt1, PetscReal *pt2) { return (T)((*this->k)(this->dim, pt1, pt2, this->ctx)); }
122f17f7ee2SSatish Balay };
123f17f7ee2SSatish Balay 
124f17f7ee2SSatish Balay   #include <../src/mat/impls/h2opus/math2opussampler.hpp>
125f17f7ee2SSatish Balay 
126f17f7ee2SSatish Balay   /* just to not clutter the code */
127f17f7ee2SSatish Balay   #if !defined(H2OPUS_USE_GPU)
128f17f7ee2SSatish Balay typedef HMatrix HMatrix_GPU;
129f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
130f17f7ee2SSatish Balay typedef DistributedHMatrix DistributedHMatrix_GPU;
131f17f7ee2SSatish Balay     #endif
132f17f7ee2SSatish Balay   #endif
133f17f7ee2SSatish Balay 
134f17f7ee2SSatish Balay typedef struct {
135f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
136f17f7ee2SSatish Balay   distributedH2OpusHandle_t handle;
137f17f7ee2SSatish Balay   #else
138f17f7ee2SSatish Balay   h2opusHandle_t handle;
139f17f7ee2SSatish Balay   #endif
140f17f7ee2SSatish Balay   /* Sequential and parallel matrices are two different classes at the moment */
141f17f7ee2SSatish Balay   HMatrix *hmatrix;
142f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
143f17f7ee2SSatish Balay   DistributedHMatrix *dist_hmatrix;
144f17f7ee2SSatish Balay   #else
145f17f7ee2SSatish Balay   HMatrix *dist_hmatrix; /* just to not clutter the code */
146f17f7ee2SSatish Balay   #endif
147f17f7ee2SSatish Balay   /* May use permutations */
148f17f7ee2SSatish Balay   PetscSF                           sf;
149f17f7ee2SSatish Balay   PetscLayout                       h2opus_rmap, h2opus_cmap;
150f17f7ee2SSatish Balay   IS                                h2opus_indexmap;
151f17f7ee2SSatish Balay   thrust::host_vector<PetscScalar> *xx, *yy;
152f17f7ee2SSatish Balay   PetscInt                          xxs, yys;
153f17f7ee2SSatish Balay   PetscBool                         multsetup;
154f17f7ee2SSatish Balay 
155f17f7ee2SSatish Balay   /* GPU */
156f17f7ee2SSatish Balay   HMatrix_GPU *hmatrix_gpu;
157f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
158f17f7ee2SSatish Balay   DistributedHMatrix_GPU *dist_hmatrix_gpu;
159f17f7ee2SSatish Balay   #else
160f17f7ee2SSatish Balay   HMatrix_GPU *dist_hmatrix_gpu; /* just to not clutter the code */
161f17f7ee2SSatish Balay   #endif
162f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
163f17f7ee2SSatish Balay   thrust::device_vector<PetscScalar> *xx_gpu, *yy_gpu;
164f17f7ee2SSatish Balay   PetscInt                            xxs_gpu, yys_gpu;
165f17f7ee2SSatish Balay   #endif
166f17f7ee2SSatish Balay 
167f17f7ee2SSatish Balay   /* construction from matvecs */
168f17f7ee2SSatish Balay   PetscMatrixSampler *sampler;
169f17f7ee2SSatish Balay   PetscBool           nativemult;
170f17f7ee2SSatish Balay 
171f17f7ee2SSatish Balay   /* Admissibility */
172f17f7ee2SSatish Balay   PetscReal eta;
173f17f7ee2SSatish Balay   PetscInt  leafsize;
174f17f7ee2SSatish Balay 
175f17f7ee2SSatish Balay   /* for dof reordering */
176f17f7ee2SSatish Balay   PetscPointCloud<PetscReal> *ptcloud;
177f17f7ee2SSatish Balay 
178f17f7ee2SSatish Balay   /* kernel for generating matrix entries */
179f17f7ee2SSatish Balay   PetscFunctionGenerator<PetscScalar> *kernel;
180f17f7ee2SSatish Balay 
181f17f7ee2SSatish Balay   /* basis orthogonalized? */
182f17f7ee2SSatish Balay   PetscBool orthogonal;
183f17f7ee2SSatish Balay 
184f17f7ee2SSatish Balay   /* customization */
185f17f7ee2SSatish Balay   PetscInt  basisord;
186f17f7ee2SSatish Balay   PetscInt  max_rank;
187f17f7ee2SSatish Balay   PetscInt  bs;
188f17f7ee2SSatish Balay   PetscReal rtol;
189f17f7ee2SSatish Balay   PetscInt  norm_max_samples;
190f17f7ee2SSatish Balay   PetscBool check_construction;
191f17f7ee2SSatish Balay   PetscBool hara_verbose;
192036e4bdcSSatish Balay   PetscBool resize;
193f17f7ee2SSatish Balay 
194f17f7ee2SSatish Balay   /* keeps track of MatScale values */
195f17f7ee2SSatish Balay   PetscScalar s;
196f17f7ee2SSatish Balay } Mat_H2OPUS;
197f17f7ee2SSatish Balay 
198d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_H2OPUS(Mat A)
199d71ae5a4SJacob Faibussowitsch {
200f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
201f17f7ee2SSatish Balay 
202f17f7ee2SSatish Balay   PetscFunctionBegin;
203f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
204f17f7ee2SSatish Balay   h2opusDestroyDistributedHandle(a->handle);
205f17f7ee2SSatish Balay   #else
206f17f7ee2SSatish Balay   h2opusDestroyHandle(a->handle);
207f17f7ee2SSatish Balay   #endif
208f17f7ee2SSatish Balay   delete a->dist_hmatrix;
209f17f7ee2SSatish Balay   delete a->hmatrix;
210f17f7ee2SSatish Balay   PetscCall(PetscSFDestroy(&a->sf));
211f17f7ee2SSatish Balay   PetscCall(PetscLayoutDestroy(&a->h2opus_rmap));
212f17f7ee2SSatish Balay   PetscCall(PetscLayoutDestroy(&a->h2opus_cmap));
213f17f7ee2SSatish Balay   PetscCall(ISDestroy(&a->h2opus_indexmap));
214f17f7ee2SSatish Balay   delete a->xx;
215f17f7ee2SSatish Balay   delete a->yy;
216f17f7ee2SSatish Balay   delete a->hmatrix_gpu;
217f17f7ee2SSatish Balay   delete a->dist_hmatrix_gpu;
218f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
219f17f7ee2SSatish Balay   delete a->xx_gpu;
220f17f7ee2SSatish Balay   delete a->yy_gpu;
221f17f7ee2SSatish Balay   #endif
222f17f7ee2SSatish Balay   delete a->sampler;
223f17f7ee2SSatish Balay   delete a->ptcloud;
224f17f7ee2SSatish Balay   delete a->kernel;
225f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdense_C", NULL));
226f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdensecuda_C", NULL));
227f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidense_C", NULL));
228f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidensecuda_C", NULL));
229f17f7ee2SSatish Balay   PetscCall(PetscObjectChangeTypeName((PetscObject)A, NULL));
230f17f7ee2SSatish Balay   PetscCall(PetscFree(A->data));
2313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
232f17f7ee2SSatish Balay }
233f17f7ee2SSatish Balay 
234d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusSetNativeMult(Mat A, PetscBool nm)
235d71ae5a4SJacob Faibussowitsch {
236f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
237f17f7ee2SSatish Balay   PetscBool   ish2opus;
238f17f7ee2SSatish Balay 
239f17f7ee2SSatish Balay   PetscFunctionBegin;
240f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
241f17f7ee2SSatish Balay   PetscValidLogicalCollectiveBool(A, nm, 2);
242f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
243f17f7ee2SSatish Balay   if (ish2opus) {
244f17f7ee2SSatish Balay     if (a->h2opus_rmap) { /* need to swap layouts for vector creation */
245f17f7ee2SSatish Balay       if ((!a->nativemult && nm) || (a->nativemult && !nm)) {
246f17f7ee2SSatish Balay         PetscLayout t;
247f17f7ee2SSatish Balay         t              = A->rmap;
248f17f7ee2SSatish Balay         A->rmap        = a->h2opus_rmap;
249f17f7ee2SSatish Balay         a->h2opus_rmap = t;
250f17f7ee2SSatish Balay         t              = A->cmap;
251f17f7ee2SSatish Balay         A->cmap        = a->h2opus_cmap;
252f17f7ee2SSatish Balay         a->h2opus_cmap = t;
253f17f7ee2SSatish Balay       }
254f17f7ee2SSatish Balay     }
255f17f7ee2SSatish Balay     a->nativemult = nm;
256f17f7ee2SSatish Balay   }
2573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
258f17f7ee2SSatish Balay }
259f17f7ee2SSatish Balay 
260d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusGetNativeMult(Mat A, PetscBool *nm)
261d71ae5a4SJacob Faibussowitsch {
262f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
263f17f7ee2SSatish Balay   PetscBool   ish2opus;
264f17f7ee2SSatish Balay 
265f17f7ee2SSatish Balay   PetscFunctionBegin;
266f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
2674f572ea9SToby Isaac   PetscAssertPointer(nm, 2);
268f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
269f17f7ee2SSatish Balay   PetscCheck(ish2opus, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name);
270f17f7ee2SSatish Balay   *nm = a->nativemult;
2713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272f17f7ee2SSatish Balay }
273f17f7ee2SSatish Balay 
274d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatNorm_H2OPUS(Mat A, NormType normtype, PetscReal *n)
275d71ae5a4SJacob Faibussowitsch {
276f17f7ee2SSatish Balay   PetscBool   ish2opus;
277f17f7ee2SSatish Balay   PetscInt    nmax = PETSC_DECIDE;
278f17f7ee2SSatish Balay   Mat_H2OPUS *a    = NULL;
279f17f7ee2SSatish Balay   PetscBool   mult = PETSC_FALSE;
280f17f7ee2SSatish Balay 
281f17f7ee2SSatish Balay   PetscFunctionBegin;
282f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
283f17f7ee2SSatish Balay   if (ish2opus) { /* set userdefine number of samples and fastpath for mult (norms are order independent) */
284f17f7ee2SSatish Balay     a = (Mat_H2OPUS *)A->data;
285f17f7ee2SSatish Balay 
286f17f7ee2SSatish Balay     nmax = a->norm_max_samples;
287f17f7ee2SSatish Balay     mult = a->nativemult;
288f17f7ee2SSatish Balay     PetscCall(MatH2OpusSetNativeMult(A, PETSC_TRUE));
289f17f7ee2SSatish Balay   } else {
290f17f7ee2SSatish Balay     PetscCall(PetscOptionsGetInt(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_approximate_norm_samples", &nmax, NULL));
291f17f7ee2SSatish Balay   }
292f17f7ee2SSatish Balay   PetscCall(MatApproximateNorm_Private(A, normtype, nmax, n));
293f17f7ee2SSatish Balay   if (a) PetscCall(MatH2OpusSetNativeMult(A, mult));
2943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
295f17f7ee2SSatish Balay }
296f17f7ee2SSatish Balay 
297d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatH2OpusResizeBuffers_Private(Mat A, PetscInt xN, PetscInt yN)
298d71ae5a4SJacob Faibussowitsch {
299f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
300f17f7ee2SSatish Balay   PetscInt    n;
301f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
302f17f7ee2SSatish Balay 
303f17f7ee2SSatish Balay   PetscFunctionBegin;
304f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
305f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
306f17f7ee2SSatish Balay   #endif
307f17f7ee2SSatish Balay   PetscCall(PetscSFGetGraph(h2opus->sf, NULL, &n, NULL, NULL));
308f17f7ee2SSatish Balay   if (boundtocpu) {
3099371c9d4SSatish Balay     if (h2opus->xxs < xN) {
3109371c9d4SSatish Balay       h2opus->xx->resize(n * xN);
3119371c9d4SSatish Balay       h2opus->xxs = xN;
3129371c9d4SSatish Balay     }
3139371c9d4SSatish Balay     if (h2opus->yys < yN) {
3149371c9d4SSatish Balay       h2opus->yy->resize(n * yN);
3159371c9d4SSatish Balay       h2opus->yys = yN;
3169371c9d4SSatish Balay     }
317f17f7ee2SSatish Balay   }
318f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
319f17f7ee2SSatish Balay   if (!boundtocpu) {
3209371c9d4SSatish Balay     if (h2opus->xxs_gpu < xN) {
3219371c9d4SSatish Balay       h2opus->xx_gpu->resize(n * xN);
3229371c9d4SSatish Balay       h2opus->xxs_gpu = xN;
3239371c9d4SSatish Balay     }
3249371c9d4SSatish Balay     if (h2opus->yys_gpu < yN) {
3259371c9d4SSatish Balay       h2opus->yy_gpu->resize(n * yN);
3269371c9d4SSatish Balay       h2opus->yys_gpu = yN;
3279371c9d4SSatish Balay     }
328f17f7ee2SSatish Balay   }
329f17f7ee2SSatish Balay   #endif
3303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
331f17f7ee2SSatish Balay }
332f17f7ee2SSatish Balay 
333d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultNKernel_H2OPUS(Mat A, PetscBool transA, Mat B, Mat C)
334d71ae5a4SJacob Faibussowitsch {
335f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
336f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
337f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle->handle;
338f17f7ee2SSatish Balay   #else
339f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle;
340f17f7ee2SSatish Balay   #endif
341f17f7ee2SSatish Balay   PetscBool    boundtocpu = PETSC_TRUE;
342f17f7ee2SSatish Balay   PetscScalar *xx, *yy, *uxx, *uyy;
343f17f7ee2SSatish Balay   PetscInt     blda, clda;
344f17f7ee2SSatish Balay   PetscMPIInt  size;
345f17f7ee2SSatish Balay   PetscSF      bsf, csf;
346f17f7ee2SSatish Balay   PetscBool    usesf = (PetscBool)(h2opus->sf && !h2opus->nativemult);
347f17f7ee2SSatish Balay 
348f17f7ee2SSatish Balay   PetscFunctionBegin;
349f17f7ee2SSatish Balay   HLibProfile::clear();
350f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
351f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
352f17f7ee2SSatish Balay   #endif
353f17f7ee2SSatish Balay   PetscCall(MatDenseGetLDA(B, &blda));
354f17f7ee2SSatish Balay   PetscCall(MatDenseGetLDA(C, &clda));
355f17f7ee2SSatish Balay   if (usesf) {
356f17f7ee2SSatish Balay     PetscInt n;
357f17f7ee2SSatish Balay 
3580dd791a8SStefano Zampini     PetscCall(MatDenseGetH2OpusStridedSF(B, h2opus->sf, &bsf));
3590dd791a8SStefano Zampini     PetscCall(MatDenseGetH2OpusStridedSF(C, h2opus->sf, &csf));
360f17f7ee2SSatish Balay 
361f17f7ee2SSatish Balay     PetscCall(MatH2OpusResizeBuffers_Private(A, B->cmap->N, C->cmap->N));
362f17f7ee2SSatish Balay     PetscCall(PetscSFGetGraph(h2opus->sf, NULL, &n, NULL, NULL));
363f17f7ee2SSatish Balay     blda = n;
364f17f7ee2SSatish Balay     clda = n;
365f17f7ee2SSatish Balay   }
366f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
367f17f7ee2SSatish Balay   if (boundtocpu) {
368f17f7ee2SSatish Balay     PetscCall(MatDenseGetArrayRead(B, (const PetscScalar **)&xx));
369f17f7ee2SSatish Balay     PetscCall(MatDenseGetArrayWrite(C, &yy));
370f17f7ee2SSatish Balay     if (usesf) {
371f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx);
372f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy);
373f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
374f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
375f17f7ee2SSatish Balay     } else {
376f17f7ee2SSatish Balay       uxx = xx;
377f17f7ee2SSatish Balay       uyy = yy;
378f17f7ee2SSatish Balay     }
379f17f7ee2SSatish Balay     if (size > 1) {
380f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
381036e4bdcSSatish Balay       PetscCheck(!transA || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
382f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
383f17f7ee2SSatish Balay       distributed_hgemv(/* transA ? H2Opus_Trans : H2Opus_NoTrans, */ h2opus->s, *h2opus->dist_hmatrix, uxx, blda, 0.0, uyy, clda, B->cmap->N, h2opus->handle);
384f17f7ee2SSatish Balay   #endif
385f17f7ee2SSatish Balay     } else {
386f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
387f17f7ee2SSatish Balay       hgemv(transA ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix, uxx, blda, 0.0, uyy, clda, B->cmap->N, handle);
388f17f7ee2SSatish Balay     }
389f17f7ee2SSatish Balay     PetscCall(MatDenseRestoreArrayRead(B, (const PetscScalar **)&xx));
390f17f7ee2SSatish Balay     if (usesf) {
391f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
392f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
393f17f7ee2SSatish Balay     }
394f17f7ee2SSatish Balay     PetscCall(MatDenseRestoreArrayWrite(C, &yy));
395f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
396f17f7ee2SSatish Balay   } else {
397f17f7ee2SSatish Balay     PetscBool ciscuda, biscuda;
398f17f7ee2SSatish Balay 
399f17f7ee2SSatish Balay     /* If not of type seqdensecuda, convert on the fly (i.e. allocate GPU memory) */
400f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)B, &biscuda, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
40148a46eb9SPierre Jolivet     if (!biscuda) PetscCall(MatConvert(B, MATDENSECUDA, MAT_INPLACE_MATRIX, &B));
402f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &ciscuda, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
403f17f7ee2SSatish Balay     if (!ciscuda) {
404f17f7ee2SSatish Balay       C->assembled = PETSC_TRUE;
405f17f7ee2SSatish Balay       PetscCall(MatConvert(C, MATDENSECUDA, MAT_INPLACE_MATRIX, &C));
406f17f7ee2SSatish Balay     }
407f17f7ee2SSatish Balay     PetscCall(MatDenseCUDAGetArrayRead(B, (const PetscScalar **)&xx));
408f17f7ee2SSatish Balay     PetscCall(MatDenseCUDAGetArrayWrite(C, &yy));
409f17f7ee2SSatish Balay     if (usesf) {
410f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx_gpu);
411f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy_gpu);
412f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
413f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(bsf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
414f17f7ee2SSatish Balay     } else {
415f17f7ee2SSatish Balay       uxx = xx;
416f17f7ee2SSatish Balay       uyy = yy;
417f17f7ee2SSatish Balay     }
418f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeBegin());
419f17f7ee2SSatish Balay     if (size > 1) {
420f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed GPU matrix");
421036e4bdcSSatish Balay       PetscCheck(!transA || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
422f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
423f17f7ee2SSatish 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);
424f17f7ee2SSatish Balay     #endif
425f17f7ee2SSatish Balay     } else {
426f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
427f17f7ee2SSatish Balay       hgemv(transA ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix_gpu, uxx, blda, 0.0, uyy, clda, B->cmap->N, handle);
428f17f7ee2SSatish Balay     }
429f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeEnd());
430f17f7ee2SSatish Balay     PetscCall(MatDenseCUDARestoreArrayRead(B, (const PetscScalar **)&xx));
431f17f7ee2SSatish Balay     if (usesf) {
432f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
433f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(csf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
434f17f7ee2SSatish Balay     }
435f17f7ee2SSatish Balay     PetscCall(MatDenseCUDARestoreArrayWrite(C, &yy));
43648a46eb9SPierre Jolivet     if (!biscuda) PetscCall(MatConvert(B, MATDENSE, MAT_INPLACE_MATRIX, &B));
43748a46eb9SPierre Jolivet     if (!ciscuda) PetscCall(MatConvert(C, MATDENSE, MAT_INPLACE_MATRIX, &C));
438f17f7ee2SSatish Balay   #endif
439f17f7ee2SSatish Balay   }
440f17f7ee2SSatish Balay   { /* log flops */
441f17f7ee2SSatish Balay     double gops, time, perf, dev;
442f17f7ee2SSatish Balay     HLibProfile::getHgemvPerf(gops, time, perf, dev);
443f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
444f17f7ee2SSatish Balay     if (boundtocpu) {
445f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
446f17f7ee2SSatish Balay     } else {
447f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
448f17f7ee2SSatish Balay     }
449f17f7ee2SSatish Balay   #else
450f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
451f17f7ee2SSatish Balay   #endif
452f17f7ee2SSatish Balay   }
4533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
454f17f7ee2SSatish Balay }
455f17f7ee2SSatish Balay 
456d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductNumeric_H2OPUS(Mat C)
457d71ae5a4SJacob Faibussowitsch {
458f17f7ee2SSatish Balay   Mat_Product *product = C->product;
459f17f7ee2SSatish Balay 
460f17f7ee2SSatish Balay   PetscFunctionBegin;
461f17f7ee2SSatish Balay   MatCheckProduct(C, 1);
462f17f7ee2SSatish Balay   switch (product->type) {
463d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AB:
464d71ae5a4SJacob Faibussowitsch     PetscCall(MatMultNKernel_H2OPUS(product->A, PETSC_FALSE, product->B, C));
465d71ae5a4SJacob Faibussowitsch     break;
466d71ae5a4SJacob Faibussowitsch   case MATPRODUCT_AtB:
467d71ae5a4SJacob Faibussowitsch     PetscCall(MatMultNKernel_H2OPUS(product->A, PETSC_TRUE, product->B, C));
468d71ae5a4SJacob Faibussowitsch     break;
469d71ae5a4SJacob Faibussowitsch   default:
470d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProduct type %s is not supported", MatProductTypes[product->type]);
471f17f7ee2SSatish Balay   }
4723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
473f17f7ee2SSatish Balay }
474f17f7ee2SSatish Balay 
475d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSymbolic_H2OPUS(Mat C)
476d71ae5a4SJacob Faibussowitsch {
477f17f7ee2SSatish Balay   Mat_Product *product = C->product;
478f17f7ee2SSatish Balay   PetscBool    cisdense;
479f17f7ee2SSatish Balay   Mat          A, B;
480f17f7ee2SSatish Balay 
481f17f7ee2SSatish Balay   PetscFunctionBegin;
482f17f7ee2SSatish Balay   MatCheckProduct(C, 1);
483f17f7ee2SSatish Balay   A = product->A;
484f17f7ee2SSatish Balay   B = product->B;
485f17f7ee2SSatish Balay   switch (product->type) {
486f17f7ee2SSatish Balay   case MATPRODUCT_AB:
487f17f7ee2SSatish Balay     PetscCall(MatSetSizes(C, A->rmap->n, B->cmap->n, A->rmap->N, B->cmap->N));
488f17f7ee2SSatish Balay     PetscCall(MatSetBlockSizesFromMats(C, product->A, product->B));
489f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATMPIDENSE, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
490f17f7ee2SSatish Balay     if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)product->B)->type_name));
491f17f7ee2SSatish Balay     PetscCall(MatSetUp(C));
492f17f7ee2SSatish Balay     break;
493f17f7ee2SSatish Balay   case MATPRODUCT_AtB:
494f17f7ee2SSatish Balay     PetscCall(MatSetSizes(C, A->cmap->n, B->cmap->n, A->cmap->N, B->cmap->N));
495f17f7ee2SSatish Balay     PetscCall(MatSetBlockSizesFromMats(C, product->A, product->B));
496f17f7ee2SSatish Balay     PetscCall(PetscObjectTypeCompareAny((PetscObject)C, &cisdense, MATSEQDENSE, MATMPIDENSE, MATSEQDENSECUDA, MATMPIDENSECUDA, ""));
497f17f7ee2SSatish Balay     if (!cisdense) PetscCall(MatSetType(C, ((PetscObject)product->B)->type_name));
498f17f7ee2SSatish Balay     PetscCall(MatSetUp(C));
499f17f7ee2SSatish Balay     break;
500d71ae5a4SJacob Faibussowitsch   default:
501d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatProduct type %s is not supported", MatProductTypes[product->type]);
502f17f7ee2SSatish Balay   }
503f17f7ee2SSatish Balay   C->ops->productsymbolic = NULL;
504f17f7ee2SSatish Balay   C->ops->productnumeric  = MatProductNumeric_H2OPUS;
5053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
506f17f7ee2SSatish Balay }
507f17f7ee2SSatish Balay 
508d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatProductSetFromOptions_H2OPUS(Mat C)
509d71ae5a4SJacob Faibussowitsch {
510f17f7ee2SSatish Balay   PetscFunctionBegin;
511f17f7ee2SSatish Balay   MatCheckProduct(C, 1);
512ad540459SPierre Jolivet   if (C->product->type == MATPRODUCT_AB || C->product->type == MATPRODUCT_AtB) C->ops->productsymbolic = MatProductSymbolic_H2OPUS;
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
514f17f7ee2SSatish Balay }
515f17f7ee2SSatish Balay 
516d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultKernel_H2OPUS(Mat A, Vec x, PetscScalar sy, Vec y, PetscBool trans)
517d71ae5a4SJacob Faibussowitsch {
518f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
519f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
520f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle->handle;
521f17f7ee2SSatish Balay   #else
522f17f7ee2SSatish Balay   h2opusHandle_t handle = h2opus->handle;
523f17f7ee2SSatish Balay   #endif
524f17f7ee2SSatish Balay   PetscBool    boundtocpu = PETSC_TRUE;
525f17f7ee2SSatish Balay   PetscInt     n;
526f17f7ee2SSatish Balay   PetscScalar *xx, *yy, *uxx, *uyy;
527f17f7ee2SSatish Balay   PetscMPIInt  size;
528f17f7ee2SSatish Balay   PetscBool    usesf = (PetscBool)(h2opus->sf && !h2opus->nativemult);
529f17f7ee2SSatish Balay 
530f17f7ee2SSatish Balay   PetscFunctionBegin;
531f17f7ee2SSatish Balay   HLibProfile::clear();
532f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
533f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
534f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
535f17f7ee2SSatish Balay   #endif
536f17f7ee2SSatish Balay   if (usesf) {
537f17f7ee2SSatish Balay     PetscCall(PetscSFGetGraph(h2opus->sf, NULL, &n, NULL, NULL));
538f17f7ee2SSatish Balay   } else n = A->rmap->n;
539f17f7ee2SSatish Balay   if (boundtocpu) {
540f17f7ee2SSatish Balay     PetscCall(VecGetArrayRead(x, (const PetscScalar **)&xx));
541f17f7ee2SSatish Balay     if (sy == 0.0) {
542f17f7ee2SSatish Balay       PetscCall(VecGetArrayWrite(y, &yy));
543f17f7ee2SSatish Balay     } else {
544f17f7ee2SSatish Balay       PetscCall(VecGetArray(y, &yy));
545f17f7ee2SSatish Balay     }
546f17f7ee2SSatish Balay     if (usesf) {
547f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx);
548f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy);
549f17f7ee2SSatish Balay 
550f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
551f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
552f17f7ee2SSatish Balay       if (sy != 0.0) {
553f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
554f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
555f17f7ee2SSatish Balay       }
556f17f7ee2SSatish Balay     } else {
557f17f7ee2SSatish Balay       uxx = xx;
558f17f7ee2SSatish Balay       uyy = yy;
559f17f7ee2SSatish Balay     }
560f17f7ee2SSatish Balay     if (size > 1) {
561f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
562036e4bdcSSatish Balay       PetscCheck(!trans || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
563f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
564f17f7ee2SSatish Balay       distributed_hgemv(/*trans ? H2Opus_Trans : H2Opus_NoTrans, */ h2opus->s, *h2opus->dist_hmatrix, uxx, n, sy, uyy, n, 1, h2opus->handle);
565f17f7ee2SSatish Balay   #endif
566f17f7ee2SSatish Balay     } else {
567f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
568f17f7ee2SSatish Balay       hgemv(trans ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix, uxx, n, sy, uyy, n, 1, handle);
569f17f7ee2SSatish Balay     }
570f17f7ee2SSatish Balay     PetscCall(VecRestoreArrayRead(x, (const PetscScalar **)&xx));
571f17f7ee2SSatish Balay     if (usesf) {
572f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
573f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
574f17f7ee2SSatish Balay     }
575f17f7ee2SSatish Balay     if (sy == 0.0) {
576f17f7ee2SSatish Balay       PetscCall(VecRestoreArrayWrite(y, &yy));
577f17f7ee2SSatish Balay     } else {
578f17f7ee2SSatish Balay       PetscCall(VecRestoreArray(y, &yy));
579f17f7ee2SSatish Balay     }
580f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
581f17f7ee2SSatish Balay   } else {
582f17f7ee2SSatish Balay     PetscCall(VecCUDAGetArrayRead(x, (const PetscScalar **)&xx));
583f17f7ee2SSatish Balay     if (sy == 0.0) {
584f17f7ee2SSatish Balay       PetscCall(VecCUDAGetArrayWrite(y, &yy));
585f17f7ee2SSatish Balay     } else {
586f17f7ee2SSatish Balay       PetscCall(VecCUDAGetArray(y, &yy));
587f17f7ee2SSatish Balay     }
588f17f7ee2SSatish Balay     if (usesf) {
589f17f7ee2SSatish Balay       uxx = MatH2OpusGetThrustPointer(*h2opus->xx_gpu);
590f17f7ee2SSatish Balay       uyy = MatH2OpusGetThrustPointer(*h2opus->yy_gpu);
591f17f7ee2SSatish Balay 
592f17f7ee2SSatish Balay       PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
593f17f7ee2SSatish Balay       PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, xx, uxx, MPI_REPLACE));
594f17f7ee2SSatish Balay       if (sy != 0.0) {
595f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
596f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(h2opus->sf, MPIU_SCALAR, yy, uyy, MPI_REPLACE));
597f17f7ee2SSatish Balay       }
598f17f7ee2SSatish Balay     } else {
599f17f7ee2SSatish Balay       uxx = xx;
600f17f7ee2SSatish Balay       uyy = yy;
601f17f7ee2SSatish Balay     }
602f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeBegin());
603f17f7ee2SSatish Balay     if (size > 1) {
604f17f7ee2SSatish Balay       PetscCheck(h2opus->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed GPU matrix");
605036e4bdcSSatish Balay       PetscCheck(!trans || A->symmetric, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "MatMultTranspose not yet coded in parallel");
606f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
607f17f7ee2SSatish Balay       distributed_hgemv(/*trans ? H2Opus_Trans : H2Opus_NoTrans, */ h2opus->s, *h2opus->dist_hmatrix_gpu, uxx, n, sy, uyy, n, 1, h2opus->handle);
608f17f7ee2SSatish Balay     #endif
609f17f7ee2SSatish Balay     } else {
610f17f7ee2SSatish Balay       PetscCheck(h2opus->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
611f17f7ee2SSatish Balay       hgemv(trans ? H2Opus_Trans : H2Opus_NoTrans, h2opus->s, *h2opus->hmatrix_gpu, uxx, n, sy, uyy, n, 1, handle);
612f17f7ee2SSatish Balay     }
613f17f7ee2SSatish Balay     PetscCall(PetscLogGpuTimeEnd());
614f17f7ee2SSatish Balay     PetscCall(VecCUDARestoreArrayRead(x, (const PetscScalar **)&xx));
615f17f7ee2SSatish Balay     if (usesf) {
616f17f7ee2SSatish Balay       PetscCall(PetscSFReduceBegin(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
617f17f7ee2SSatish Balay       PetscCall(PetscSFReduceEnd(h2opus->sf, MPIU_SCALAR, uyy, yy, MPI_REPLACE));
618f17f7ee2SSatish Balay     }
619f17f7ee2SSatish Balay     if (sy == 0.0) {
620f17f7ee2SSatish Balay       PetscCall(VecCUDARestoreArrayWrite(y, &yy));
621f17f7ee2SSatish Balay     } else {
622f17f7ee2SSatish Balay       PetscCall(VecCUDARestoreArray(y, &yy));
623f17f7ee2SSatish Balay     }
624f17f7ee2SSatish Balay   #endif
625f17f7ee2SSatish Balay   }
626f17f7ee2SSatish Balay   { /* log flops */
627f17f7ee2SSatish Balay     double gops, time, perf, dev;
628f17f7ee2SSatish Balay     HLibProfile::getHgemvPerf(gops, time, perf, dev);
629f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
630f17f7ee2SSatish Balay     if (boundtocpu) {
631f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
632f17f7ee2SSatish Balay     } else {
633f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
634f17f7ee2SSatish Balay     }
635f17f7ee2SSatish Balay   #else
636f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
637f17f7ee2SSatish Balay   #endif
638f17f7ee2SSatish Balay   }
6393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
640f17f7ee2SSatish Balay }
641f17f7ee2SSatish Balay 
642d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_H2OPUS(Mat A, Vec x, Vec y)
643d71ae5a4SJacob Faibussowitsch {
644f17f7ee2SSatish Balay   PetscBool xiscuda, yiscuda;
645f17f7ee2SSatish Balay 
646f17f7ee2SSatish Balay   PetscFunctionBegin;
647f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
648f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)y, &yiscuda, VECSEQCUDA, VECMPICUDA, ""));
649f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !yiscuda));
650f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 0.0, y, PETSC_TRUE));
6513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
652f17f7ee2SSatish Balay }
653f17f7ee2SSatish Balay 
654d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_H2OPUS(Mat A, Vec x, Vec y)
655d71ae5a4SJacob Faibussowitsch {
656f17f7ee2SSatish Balay   PetscBool xiscuda, yiscuda;
657f17f7ee2SSatish Balay 
658f17f7ee2SSatish Balay   PetscFunctionBegin;
659f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
660f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)y, &yiscuda, VECSEQCUDA, VECMPICUDA, ""));
661f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !yiscuda));
662f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 0.0, y, PETSC_FALSE));
6633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
664f17f7ee2SSatish Balay }
665f17f7ee2SSatish Balay 
666d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_H2OPUS(Mat A, Vec x, Vec y, Vec z)
667d71ae5a4SJacob Faibussowitsch {
668f17f7ee2SSatish Balay   PetscBool xiscuda, ziscuda;
669f17f7ee2SSatish Balay 
670f17f7ee2SSatish Balay   PetscFunctionBegin;
671f17f7ee2SSatish Balay   PetscCall(VecCopy(y, z));
672f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
673f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)z, &ziscuda, VECSEQCUDA, VECMPICUDA, ""));
674f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !ziscuda));
675f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 1.0, z, PETSC_TRUE));
6763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
677f17f7ee2SSatish Balay }
678f17f7ee2SSatish Balay 
679d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_H2OPUS(Mat A, Vec x, Vec y, Vec z)
680d71ae5a4SJacob Faibussowitsch {
681f17f7ee2SSatish Balay   PetscBool xiscuda, ziscuda;
682f17f7ee2SSatish Balay 
683f17f7ee2SSatish Balay   PetscFunctionBegin;
684f17f7ee2SSatish Balay   PetscCall(VecCopy(y, z));
685f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)x, &xiscuda, VECSEQCUDA, VECMPICUDA, ""));
686f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompareAny((PetscObject)z, &ziscuda, VECSEQCUDA, VECMPICUDA, ""));
687f17f7ee2SSatish Balay   PetscCall(MatH2OpusUpdateIfNeeded(A, !xiscuda || !ziscuda));
688f17f7ee2SSatish Balay   PetscCall(MatMultKernel_H2OPUS(A, x, 1.0, z, PETSC_FALSE));
6893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
690f17f7ee2SSatish Balay }
691f17f7ee2SSatish Balay 
692d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_H2OPUS(Mat A, PetscScalar s)
693d71ae5a4SJacob Faibussowitsch {
694f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
695f17f7ee2SSatish Balay 
696f17f7ee2SSatish Balay   PetscFunctionBegin;
697f17f7ee2SSatish Balay   a->s *= s;
6983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
699f17f7ee2SSatish Balay }
700f17f7ee2SSatish Balay 
701d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetFromOptions_H2OPUS(Mat A, PetscOptionItems *PetscOptionsObject)
702d71ae5a4SJacob Faibussowitsch {
703f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
704f17f7ee2SSatish Balay 
705f17f7ee2SSatish Balay   PetscFunctionBegin;
706036e4bdcSSatish Balay   PetscOptionsHeadBegin(PetscOptionsObject, "H2OPUS options");
707f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_leafsize", "Leaf size of cluster tree", NULL, a->leafsize, &a->leafsize, NULL));
708f17f7ee2SSatish Balay   PetscCall(PetscOptionsReal("-mat_h2opus_eta", "Admissibility condition tolerance", NULL, a->eta, &a->eta, NULL));
709f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_order", "Basis order for off-diagonal sampling when constructed from kernel", NULL, a->basisord, &a->basisord, NULL));
710f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_maxrank", "Maximum rank when constructed from matvecs", NULL, a->max_rank, &a->max_rank, NULL));
711f17f7ee2SSatish Balay   PetscCall(PetscOptionsInt("-mat_h2opus_samples", "Maximum number of samples to be taken concurrently when constructing from matvecs", NULL, a->bs, &a->bs, NULL));
712aaa8cc7dSPierre Jolivet   PetscCall(PetscOptionsInt("-mat_h2opus_normsamples", "Maximum number of samples to be when estimating norms", NULL, a->norm_max_samples, &a->norm_max_samples, NULL));
713f17f7ee2SSatish Balay   PetscCall(PetscOptionsReal("-mat_h2opus_rtol", "Relative tolerance for construction from sampling", NULL, a->rtol, &a->rtol, NULL));
714f17f7ee2SSatish Balay   PetscCall(PetscOptionsBool("-mat_h2opus_check", "Check error when constructing from sampling during MatAssemblyEnd()", NULL, a->check_construction, &a->check_construction, NULL));
715f17f7ee2SSatish Balay   PetscCall(PetscOptionsBool("-mat_h2opus_hara_verbose", "Verbose output from hara construction", NULL, a->hara_verbose, &a->hara_verbose, NULL));
716036e4bdcSSatish Balay   PetscCall(PetscOptionsBool("-mat_h2opus_resize", "Resize after compression", NULL, a->resize, &a->resize, NULL));
717036e4bdcSSatish Balay   PetscOptionsHeadEnd();
7183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
719f17f7ee2SSatish Balay }
720f17f7ee2SSatish Balay 
7218434afd1SBarry Smith static PetscErrorCode MatH2OpusSetCoords_H2OPUS(Mat, PetscInt, const PetscReal[], PetscBool, MatH2OpusKernelFn *, void *);
722f17f7ee2SSatish Balay 
723d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatH2OpusInferCoordinates_Private(Mat A)
724d71ae5a4SJacob Faibussowitsch {
725f17f7ee2SSatish Balay   Mat_H2OPUS        *a = (Mat_H2OPUS *)A->data;
726f17f7ee2SSatish Balay   Vec                c;
727f17f7ee2SSatish Balay   PetscInt           spacedim;
728f17f7ee2SSatish Balay   const PetscScalar *coords;
729f17f7ee2SSatish Balay 
730f17f7ee2SSatish Balay   PetscFunctionBegin;
7313ba16761SJacob Faibussowitsch   if (a->ptcloud) PetscFunctionReturn(PETSC_SUCCESS);
732f17f7ee2SSatish Balay   PetscCall(PetscObjectQuery((PetscObject)A, "__math2opus_coords", (PetscObject *)&c));
733f17f7ee2SSatish Balay   if (!c && a->sampler) {
734f17f7ee2SSatish Balay     Mat S = a->sampler->GetSamplingMat();
735f17f7ee2SSatish Balay 
736f17f7ee2SSatish Balay     PetscCall(PetscObjectQuery((PetscObject)S, "__math2opus_coords", (PetscObject *)&c));
737f17f7ee2SSatish Balay   }
738f17f7ee2SSatish Balay   if (!c) {
739f17f7ee2SSatish Balay     PetscCall(MatH2OpusSetCoords_H2OPUS(A, -1, NULL, PETSC_FALSE, NULL, NULL));
740f17f7ee2SSatish Balay   } else {
741f17f7ee2SSatish Balay     PetscCall(VecGetArrayRead(c, &coords));
742f17f7ee2SSatish Balay     PetscCall(VecGetBlockSize(c, &spacedim));
743f17f7ee2SSatish Balay     PetscCall(MatH2OpusSetCoords_H2OPUS(A, spacedim, coords, PETSC_FALSE, NULL, NULL));
744f17f7ee2SSatish Balay     PetscCall(VecRestoreArrayRead(c, &coords));
745f17f7ee2SSatish Balay   }
7463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
747f17f7ee2SSatish Balay }
748f17f7ee2SSatish Balay 
749d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSetUpMultiply_H2OPUS(Mat A)
750d71ae5a4SJacob Faibussowitsch {
751f17f7ee2SSatish Balay   MPI_Comm      comm;
752f17f7ee2SSatish Balay   PetscMPIInt   size;
753f17f7ee2SSatish Balay   Mat_H2OPUS   *a = (Mat_H2OPUS *)A->data;
754f17f7ee2SSatish Balay   PetscInt      n = 0, *idx = NULL;
755f17f7ee2SSatish Balay   int          *iidx = NULL;
756f17f7ee2SSatish Balay   PetscCopyMode own;
757f17f7ee2SSatish Balay   PetscBool     rid;
758f17f7ee2SSatish Balay 
759f17f7ee2SSatish Balay   PetscFunctionBegin;
7603ba16761SJacob Faibussowitsch   if (a->multsetup) PetscFunctionReturn(PETSC_SUCCESS);
761f17f7ee2SSatish Balay   if (a->sf) { /* MatDuplicate_H2OPUS takes reference to the SF */
762f17f7ee2SSatish Balay     PetscCall(PetscSFGetGraph(a->sf, NULL, &n, NULL, NULL));
763f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
764f17f7ee2SSatish Balay     a->xx_gpu  = new thrust::device_vector<PetscScalar>(n);
765f17f7ee2SSatish Balay     a->yy_gpu  = new thrust::device_vector<PetscScalar>(n);
766f17f7ee2SSatish Balay     a->xxs_gpu = 1;
767f17f7ee2SSatish Balay     a->yys_gpu = 1;
768f17f7ee2SSatish Balay   #endif
769f17f7ee2SSatish Balay     a->xx  = new thrust::host_vector<PetscScalar>(n);
770f17f7ee2SSatish Balay     a->yy  = new thrust::host_vector<PetscScalar>(n);
771f17f7ee2SSatish Balay     a->xxs = 1;
772f17f7ee2SSatish Balay     a->yys = 1;
773f17f7ee2SSatish Balay   } else {
774f17f7ee2SSatish Balay     IS is;
775f17f7ee2SSatish Balay     PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
776f17f7ee2SSatish Balay     PetscCallMPI(MPI_Comm_size(comm, &size));
777f17f7ee2SSatish Balay     if (!a->h2opus_indexmap) {
778f17f7ee2SSatish Balay       if (size > 1) {
779f17f7ee2SSatish Balay         PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
780f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
781f17f7ee2SSatish Balay         iidx = MatH2OpusGetThrustPointer(a->dist_hmatrix->basis_tree.basis_branch.index_map);
782f17f7ee2SSatish Balay         n    = a->dist_hmatrix->basis_tree.basis_branch.index_map.size();
783f17f7ee2SSatish Balay   #endif
784f17f7ee2SSatish Balay       } else {
785f17f7ee2SSatish Balay         iidx = MatH2OpusGetThrustPointer(a->hmatrix->u_basis_tree.index_map);
786f17f7ee2SSatish Balay         n    = a->hmatrix->u_basis_tree.index_map.size();
787f17f7ee2SSatish Balay       }
788f17f7ee2SSatish Balay 
789f17f7ee2SSatish Balay       if (PetscDefined(USE_64BIT_INDICES)) {
790f17f7ee2SSatish Balay         PetscInt i;
791f17f7ee2SSatish Balay 
792f17f7ee2SSatish Balay         own = PETSC_OWN_POINTER;
793f17f7ee2SSatish Balay         PetscCall(PetscMalloc1(n, &idx));
794f17f7ee2SSatish Balay         for (i = 0; i < n; i++) idx[i] = iidx[i];
795f17f7ee2SSatish Balay       } else {
796f17f7ee2SSatish Balay         own = PETSC_COPY_VALUES;
797f17f7ee2SSatish Balay         idx = (PetscInt *)iidx;
798f17f7ee2SSatish Balay       }
799f17f7ee2SSatish Balay       PetscCall(ISCreateGeneral(comm, n, idx, own, &is));
800f17f7ee2SSatish Balay       PetscCall(ISSetPermutation(is));
801f17f7ee2SSatish Balay       PetscCall(ISViewFromOptions(is, (PetscObject)A, "-mat_h2opus_indexmap_view"));
802f17f7ee2SSatish Balay       a->h2opus_indexmap = is;
803f17f7ee2SSatish Balay     }
804f17f7ee2SSatish Balay     PetscCall(ISGetLocalSize(a->h2opus_indexmap, &n));
805f17f7ee2SSatish Balay     PetscCall(ISGetIndices(a->h2opus_indexmap, (const PetscInt **)&idx));
806f17f7ee2SSatish Balay     rid = (PetscBool)(n == A->rmap->n);
807f17f7ee2SSatish Balay     PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &rid, 1, MPIU_BOOL, MPI_LAND, comm));
80848a46eb9SPierre Jolivet     if (rid) PetscCall(ISIdentity(a->h2opus_indexmap, &rid));
809f17f7ee2SSatish Balay     if (!rid) {
810f17f7ee2SSatish Balay       if (size > 1) { /* Parallel distribution may be different, save it here for fast path in MatMult (see MatH2OpusSetNativeMult) */
811f17f7ee2SSatish Balay         PetscCall(PetscLayoutCreate(comm, &a->h2opus_rmap));
812f17f7ee2SSatish Balay         PetscCall(PetscLayoutSetLocalSize(a->h2opus_rmap, n));
813f17f7ee2SSatish Balay         PetscCall(PetscLayoutSetUp(a->h2opus_rmap));
814f17f7ee2SSatish Balay         PetscCall(PetscLayoutReference(a->h2opus_rmap, &a->h2opus_cmap));
815f17f7ee2SSatish Balay       }
816f17f7ee2SSatish Balay       PetscCall(PetscSFCreate(comm, &a->sf));
817f17f7ee2SSatish Balay       PetscCall(PetscSFSetGraphLayout(a->sf, A->rmap, n, NULL, PETSC_OWN_POINTER, idx));
818036e4bdcSSatish Balay       PetscCall(PetscSFSetUp(a->sf));
819f17f7ee2SSatish Balay       PetscCall(PetscSFViewFromOptions(a->sf, (PetscObject)A, "-mat_h2opus_sf_view"));
820f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
821f17f7ee2SSatish Balay       a->xx_gpu  = new thrust::device_vector<PetscScalar>(n);
822f17f7ee2SSatish Balay       a->yy_gpu  = new thrust::device_vector<PetscScalar>(n);
823f17f7ee2SSatish Balay       a->xxs_gpu = 1;
824f17f7ee2SSatish Balay       a->yys_gpu = 1;
825f17f7ee2SSatish Balay   #endif
826f17f7ee2SSatish Balay       a->xx  = new thrust::host_vector<PetscScalar>(n);
827f17f7ee2SSatish Balay       a->yy  = new thrust::host_vector<PetscScalar>(n);
828f17f7ee2SSatish Balay       a->xxs = 1;
829f17f7ee2SSatish Balay       a->yys = 1;
830f17f7ee2SSatish Balay     }
831f17f7ee2SSatish Balay     PetscCall(ISRestoreIndices(a->h2opus_indexmap, (const PetscInt **)&idx));
832f17f7ee2SSatish Balay   }
833f17f7ee2SSatish Balay   a->multsetup = PETSC_TRUE;
8343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
835f17f7ee2SSatish Balay }
836f17f7ee2SSatish Balay 
837d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_H2OPUS(Mat A, MatAssemblyType assemblytype)
838d71ae5a4SJacob Faibussowitsch {
839f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
840f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
841f17f7ee2SSatish Balay   h2opusHandle_t handle = a->handle->handle;
842f17f7ee2SSatish Balay   #else
843f17f7ee2SSatish Balay   h2opusHandle_t handle = a->handle;
844f17f7ee2SSatish Balay   #endif
845f17f7ee2SSatish Balay   PetscBool   kernel       = PETSC_FALSE;
846f17f7ee2SSatish Balay   PetscBool   boundtocpu   = PETSC_TRUE;
847f17f7ee2SSatish Balay   PetscBool   samplingdone = PETSC_FALSE;
848f17f7ee2SSatish Balay   MPI_Comm    comm;
849f17f7ee2SSatish Balay   PetscMPIInt size;
850f17f7ee2SSatish Balay 
851f17f7ee2SSatish Balay   PetscFunctionBegin;
852f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
853036e4bdcSSatish Balay   PetscCheck(A->rmap->n == A->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Different row and column local sizes are not supported");
854036e4bdcSSatish Balay   PetscCheck(A->rmap->N == A->cmap->N, comm, PETSC_ERR_SUP, "Rectangular matrices are not supported");
855f17f7ee2SSatish Balay 
856f17f7ee2SSatish Balay   /* XXX */
857f17f7ee2SSatish Balay   a->leafsize = PetscMin(a->leafsize, PetscMin(A->rmap->N, A->cmap->N));
858f17f7ee2SSatish Balay 
859f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(comm, &size));
860f17f7ee2SSatish Balay   /* TODO REUSABILITY of geometric construction */
861f17f7ee2SSatish Balay   delete a->hmatrix;
862f17f7ee2SSatish Balay   delete a->dist_hmatrix;
863f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
864f17f7ee2SSatish Balay   delete a->hmatrix_gpu;
865f17f7ee2SSatish Balay   delete a->dist_hmatrix_gpu;
866f17f7ee2SSatish Balay   #endif
867f17f7ee2SSatish Balay   a->orthogonal = PETSC_FALSE;
868f17f7ee2SSatish Balay 
869f17f7ee2SSatish Balay   /* TODO: other? */
870f17f7ee2SSatish Balay   H2OpusBoxCenterAdmissibility adm(a->eta);
871f17f7ee2SSatish Balay 
872f17f7ee2SSatish Balay   PetscCall(PetscLogEventBegin(MAT_H2Opus_Build, A, 0, 0, 0));
873f17f7ee2SSatish Balay   if (size > 1) {
874f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
875f17f7ee2SSatish Balay     a->dist_hmatrix = new DistributedHMatrix(A->rmap->n /* ,A->symmetric */);
876f17f7ee2SSatish Balay   #else
877f17f7ee2SSatish Balay     a->dist_hmatrix = NULL;
878f17f7ee2SSatish Balay   #endif
879b94d7dedSBarry Smith   } else a->hmatrix = new HMatrix(A->rmap->n, A->symmetric == PETSC_BOOL3_TRUE);
880f17f7ee2SSatish Balay   PetscCall(MatH2OpusInferCoordinates_Private(A));
881f17f7ee2SSatish Balay   PetscCheck(a->ptcloud, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing pointcloud");
882f17f7ee2SSatish Balay   if (a->kernel) {
883f17f7ee2SSatish Balay     BoxEntryGen<PetscScalar, H2OPUS_HWTYPE_CPU, PetscFunctionGenerator<PetscScalar>> entry_gen(*a->kernel);
884f17f7ee2SSatish Balay     if (size > 1) {
885f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
886f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
887f17f7ee2SSatish Balay       buildDistributedHMatrix(*a->dist_hmatrix, a->ptcloud, adm, entry_gen, a->leafsize, a->basisord, a->handle);
888f17f7ee2SSatish Balay   #endif
889f17f7ee2SSatish Balay     } else {
890f17f7ee2SSatish Balay       buildHMatrix(*a->hmatrix, a->ptcloud, adm, entry_gen, a->leafsize, a->basisord);
891f17f7ee2SSatish Balay     }
892f17f7ee2SSatish Balay     kernel = PETSC_TRUE;
893f17f7ee2SSatish Balay   } else {
894036e4bdcSSatish Balay     PetscCheck(size <= 1, comm, PETSC_ERR_SUP, "Construction from sampling not supported in parallel");
895f17f7ee2SSatish Balay     buildHMatrixStructure(*a->hmatrix, a->ptcloud, a->leafsize, adm);
896f17f7ee2SSatish Balay   }
897f17f7ee2SSatish Balay   PetscCall(MatSetUpMultiply_H2OPUS(A));
898f17f7ee2SSatish Balay 
899f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
900f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
901f17f7ee2SSatish Balay   if (!boundtocpu) {
902f17f7ee2SSatish Balay     if (size > 1) {
903f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing distributed CPU matrix");
904f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
905f17f7ee2SSatish Balay       a->dist_hmatrix_gpu = new DistributedHMatrix_GPU(*a->dist_hmatrix);
906f17f7ee2SSatish Balay     #endif
907f17f7ee2SSatish Balay     } else {
908f17f7ee2SSatish Balay       a->hmatrix_gpu = new HMatrix_GPU(*a->hmatrix);
909f17f7ee2SSatish Balay     }
910f17f7ee2SSatish Balay   }
911f17f7ee2SSatish Balay   #endif
912f17f7ee2SSatish Balay   if (size == 1) {
913f17f7ee2SSatish Balay     if (!kernel && a->sampler && a->sampler->GetSamplingMat()) {
914f17f7ee2SSatish Balay       PetscReal Anorm;
915f17f7ee2SSatish Balay       bool      verbose;
916f17f7ee2SSatish Balay 
917f17f7ee2SSatish Balay       PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_hara_verbose", &a->hara_verbose, NULL));
918f17f7ee2SSatish Balay       verbose = a->hara_verbose;
919f17f7ee2SSatish Balay       PetscCall(MatApproximateNorm_Private(a->sampler->GetSamplingMat(), NORM_2, a->norm_max_samples, &Anorm));
920f17f7ee2SSatish 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));
921ad540459SPierre 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());
922f17f7ee2SSatish Balay       a->sampler->SetStream(handle->getMainStream());
923f17f7ee2SSatish Balay       if (boundtocpu) {
924f17f7ee2SSatish Balay         a->sampler->SetGPUSampling(false);
925f17f7ee2SSatish Balay         hara(a->sampler, *a->hmatrix, a->max_rank, 10 /* TODO */, a->rtol * Anorm, a->bs, handle, verbose);
926f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
927f17f7ee2SSatish Balay       } else {
928f17f7ee2SSatish Balay         a->sampler->SetGPUSampling(true);
929f17f7ee2SSatish Balay         hara(a->sampler, *a->hmatrix_gpu, a->max_rank, 10 /* TODO */, a->rtol * Anorm, a->bs, handle, verbose);
930f17f7ee2SSatish Balay   #endif
931f17f7ee2SSatish Balay       }
932f17f7ee2SSatish Balay       samplingdone = PETSC_TRUE;
933f17f7ee2SSatish Balay     }
934f17f7ee2SSatish Balay   }
935f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
936f17f7ee2SSatish Balay   if (!boundtocpu) {
937f17f7ee2SSatish Balay     delete a->hmatrix;
938f17f7ee2SSatish Balay     delete a->dist_hmatrix;
939f17f7ee2SSatish Balay     a->hmatrix      = NULL;
940f17f7ee2SSatish Balay     a->dist_hmatrix = NULL;
941f17f7ee2SSatish Balay   }
942f17f7ee2SSatish Balay   A->offloadmask = boundtocpu ? PETSC_OFFLOAD_CPU : PETSC_OFFLOAD_GPU;
943f17f7ee2SSatish Balay   #endif
944f17f7ee2SSatish Balay   PetscCall(PetscLogEventEnd(MAT_H2Opus_Build, A, 0, 0, 0));
945f17f7ee2SSatish Balay 
946f17f7ee2SSatish Balay   if (!a->s) a->s = 1.0;
947f17f7ee2SSatish Balay   A->assembled = PETSC_TRUE;
948f17f7ee2SSatish Balay 
949f17f7ee2SSatish Balay   if (samplingdone) {
950f17f7ee2SSatish Balay     PetscBool check  = a->check_construction;
951f17f7ee2SSatish Balay     PetscBool checke = PETSC_FALSE;
952f17f7ee2SSatish Balay 
953f17f7ee2SSatish Balay     PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_check", &check, NULL));
954f17f7ee2SSatish Balay     PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_check_explicit", &checke, NULL));
955f17f7ee2SSatish Balay     if (check) {
956f17f7ee2SSatish Balay       Mat       E, Ae;
957f17f7ee2SSatish Balay       PetscReal n1, ni, n2;
958f17f7ee2SSatish Balay       PetscReal n1A, niA, n2A;
959f17f7ee2SSatish Balay       void (*normfunc)(void);
960f17f7ee2SSatish Balay 
961f17f7ee2SSatish Balay       Ae = a->sampler->GetSamplingMat();
962f17f7ee2SSatish Balay       PetscCall(MatConvert(A, MATSHELL, MAT_INITIAL_MATRIX, &E));
963f17f7ee2SSatish Balay       PetscCall(MatShellSetOperation(E, MATOP_NORM, (void (*)(void))MatNorm_H2OPUS));
964f17f7ee2SSatish Balay       PetscCall(MatAXPY(E, -1.0, Ae, DIFFERENT_NONZERO_PATTERN));
965f17f7ee2SSatish Balay       PetscCall(MatNorm(E, NORM_1, &n1));
966f17f7ee2SSatish Balay       PetscCall(MatNorm(E, NORM_INFINITY, &ni));
967f17f7ee2SSatish Balay       PetscCall(MatNorm(E, NORM_2, &n2));
968f17f7ee2SSatish Balay       if (checke) {
969f17f7ee2SSatish Balay         Mat eA, eE, eAe;
970f17f7ee2SSatish Balay 
971f17f7ee2SSatish Balay         PetscCall(MatComputeOperator(A, MATAIJ, &eA));
972f17f7ee2SSatish Balay         PetscCall(MatComputeOperator(E, MATAIJ, &eE));
973f17f7ee2SSatish Balay         PetscCall(MatComputeOperator(Ae, MATAIJ, &eAe));
9742ce66baaSPierre Jolivet         PetscCall(MatFilter(eA, PETSC_SMALL, PETSC_FALSE, PETSC_FALSE));
9752ce66baaSPierre Jolivet         PetscCall(MatFilter(eE, PETSC_SMALL, PETSC_FALSE, PETSC_FALSE));
9762ce66baaSPierre Jolivet         PetscCall(MatFilter(eAe, PETSC_SMALL, PETSC_FALSE, PETSC_FALSE));
977f17f7ee2SSatish Balay         PetscCall(PetscObjectSetName((PetscObject)eA, "H2Mat"));
978f17f7ee2SSatish Balay         PetscCall(MatView(eA, NULL));
979f17f7ee2SSatish Balay         PetscCall(PetscObjectSetName((PetscObject)eAe, "S"));
980f17f7ee2SSatish Balay         PetscCall(MatView(eAe, NULL));
981f17f7ee2SSatish Balay         PetscCall(PetscObjectSetName((PetscObject)eE, "H2Mat - S"));
982f17f7ee2SSatish Balay         PetscCall(MatView(eE, NULL));
983f17f7ee2SSatish Balay         PetscCall(MatDestroy(&eA));
984f17f7ee2SSatish Balay         PetscCall(MatDestroy(&eE));
985f17f7ee2SSatish Balay         PetscCall(MatDestroy(&eAe));
986f17f7ee2SSatish Balay       }
987f17f7ee2SSatish Balay 
988f17f7ee2SSatish Balay       PetscCall(MatGetOperation(Ae, MATOP_NORM, &normfunc));
989f17f7ee2SSatish Balay       PetscCall(MatSetOperation(Ae, MATOP_NORM, (void (*)(void))MatNorm_H2OPUS));
990f17f7ee2SSatish Balay       PetscCall(MatNorm(Ae, NORM_1, &n1A));
991f17f7ee2SSatish Balay       PetscCall(MatNorm(Ae, NORM_INFINITY, &niA));
992f17f7ee2SSatish Balay       PetscCall(MatNorm(Ae, NORM_2, &n2A));
993f17f7ee2SSatish Balay       n1A = PetscMax(n1A, PETSC_SMALL);
994f17f7ee2SSatish Balay       n2A = PetscMax(n2A, PETSC_SMALL);
995f17f7ee2SSatish Balay       niA = PetscMax(niA, PETSC_SMALL);
996f17f7ee2SSatish Balay       PetscCall(MatSetOperation(Ae, MATOP_NORM, normfunc));
997f17f7ee2SSatish 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)));
998f17f7ee2SSatish Balay       PetscCall(MatDestroy(&E));
999f17f7ee2SSatish Balay     }
1000f17f7ee2SSatish Balay     a->sampler->SetSamplingMat(NULL);
1001f17f7ee2SSatish Balay   }
10023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1003f17f7ee2SSatish Balay }
1004f17f7ee2SSatish Balay 
1005d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_H2OPUS(Mat A)
1006d71ae5a4SJacob Faibussowitsch {
1007f17f7ee2SSatish Balay   PetscMPIInt size;
1008f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1009f17f7ee2SSatish Balay 
1010f17f7ee2SSatish Balay   PetscFunctionBegin;
1011f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1012036e4bdcSSatish Balay   PetscCheck(size <= 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not yet supported");
1013f17f7ee2SSatish Balay   a->hmatrix->clearData();
1014f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1015f17f7ee2SSatish Balay   if (a->hmatrix_gpu) a->hmatrix_gpu->clearData();
1016f17f7ee2SSatish Balay   #endif
10173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1018f17f7ee2SSatish Balay }
1019f17f7ee2SSatish Balay 
1020d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_H2OPUS(Mat B, MatDuplicateOption op, Mat *nA)
1021d71ae5a4SJacob Faibussowitsch {
1022f17f7ee2SSatish Balay   Mat         A;
1023f17f7ee2SSatish Balay   Mat_H2OPUS *a, *b = (Mat_H2OPUS *)B->data;
1024f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1025f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_FALSE;
1026f17f7ee2SSatish Balay   #else
1027f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_TRUE;
1028f17f7ee2SSatish Balay   #endif
1029f17f7ee2SSatish Balay   MPI_Comm comm;
1030f17f7ee2SSatish Balay 
1031f17f7ee2SSatish Balay   PetscFunctionBegin;
1032f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)B, &comm));
1033f17f7ee2SSatish Balay   PetscCall(MatCreate(comm, &A));
1034f17f7ee2SSatish Balay   PetscCall(MatSetSizes(A, B->rmap->n, B->cmap->n, B->rmap->N, B->cmap->N));
1035f17f7ee2SSatish Balay   PetscCall(MatSetType(A, MATH2OPUS));
1036f17f7ee2SSatish Balay   PetscCall(MatPropagateSymmetryOptions(B, A));
1037f17f7ee2SSatish Balay   a = (Mat_H2OPUS *)A->data;
1038f17f7ee2SSatish Balay 
1039f17f7ee2SSatish Balay   a->eta              = b->eta;
1040f17f7ee2SSatish Balay   a->leafsize         = b->leafsize;
1041f17f7ee2SSatish Balay   a->basisord         = b->basisord;
1042f17f7ee2SSatish Balay   a->max_rank         = b->max_rank;
1043f17f7ee2SSatish Balay   a->bs               = b->bs;
1044f17f7ee2SSatish Balay   a->rtol             = b->rtol;
1045f17f7ee2SSatish Balay   a->norm_max_samples = b->norm_max_samples;
1046f17f7ee2SSatish Balay   if (op == MAT_COPY_VALUES) a->s = b->s;
1047f17f7ee2SSatish Balay 
1048f17f7ee2SSatish Balay   a->ptcloud = new PetscPointCloud<PetscReal>(*b->ptcloud);
1049f17f7ee2SSatish Balay   if (op == MAT_COPY_VALUES && b->kernel) a->kernel = new PetscFunctionGenerator<PetscScalar>(*b->kernel);
1050f17f7ee2SSatish Balay 
1051f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1052ad540459SPierre Jolivet   if (b->dist_hmatrix) a->dist_hmatrix = new DistributedHMatrix(*b->dist_hmatrix);
1053f17f7ee2SSatish Balay     #if defined(PETSC_H2OPUS_USE_GPU)
1054ad540459SPierre Jolivet   if (b->dist_hmatrix_gpu) a->dist_hmatrix_gpu = new DistributedHMatrix_GPU(*b->dist_hmatrix_gpu);
1055f17f7ee2SSatish Balay     #endif
1056f17f7ee2SSatish Balay   #endif
1057f17f7ee2SSatish Balay   if (b->hmatrix) {
1058f17f7ee2SSatish Balay     a->hmatrix = new HMatrix(*b->hmatrix);
1059f17f7ee2SSatish Balay     if (op == MAT_DO_NOT_COPY_VALUES) a->hmatrix->clearData();
1060f17f7ee2SSatish Balay   }
1061f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1062f17f7ee2SSatish Balay   if (b->hmatrix_gpu) {
1063f17f7ee2SSatish Balay     a->hmatrix_gpu = new HMatrix_GPU(*b->hmatrix_gpu);
1064f17f7ee2SSatish Balay     if (op == MAT_DO_NOT_COPY_VALUES) a->hmatrix_gpu->clearData();
1065f17f7ee2SSatish Balay   }
1066f17f7ee2SSatish Balay   #endif
1067f17f7ee2SSatish Balay   if (b->sf) {
1068f17f7ee2SSatish Balay     PetscCall(PetscObjectReference((PetscObject)b->sf));
1069f17f7ee2SSatish Balay     a->sf = b->sf;
1070f17f7ee2SSatish Balay   }
1071f17f7ee2SSatish Balay   if (b->h2opus_indexmap) {
1072f17f7ee2SSatish Balay     PetscCall(PetscObjectReference((PetscObject)b->h2opus_indexmap));
1073f17f7ee2SSatish Balay     a->h2opus_indexmap = b->h2opus_indexmap;
1074f17f7ee2SSatish Balay   }
1075f17f7ee2SSatish Balay 
1076f17f7ee2SSatish Balay   PetscCall(MatSetUp(A));
1077f17f7ee2SSatish Balay   PetscCall(MatSetUpMultiply_H2OPUS(A));
1078f17f7ee2SSatish Balay   if (op == MAT_COPY_VALUES) {
1079f17f7ee2SSatish Balay     A->assembled  = PETSC_TRUE;
1080f17f7ee2SSatish Balay     a->orthogonal = b->orthogonal;
1081f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1082f17f7ee2SSatish Balay     A->offloadmask = B->offloadmask;
1083f17f7ee2SSatish Balay   #endif
1084f17f7ee2SSatish Balay   }
1085f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1086f17f7ee2SSatish Balay   iscpu = B->boundtocpu;
1087f17f7ee2SSatish Balay   #endif
1088f17f7ee2SSatish Balay   PetscCall(MatBindToCPU(A, iscpu));
1089f17f7ee2SSatish Balay 
1090f17f7ee2SSatish Balay   *nA = A;
10913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1092f17f7ee2SSatish Balay }
1093f17f7ee2SSatish Balay 
1094d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_H2OPUS(Mat A, PetscViewer view)
1095d71ae5a4SJacob Faibussowitsch {
1096f17f7ee2SSatish Balay   Mat_H2OPUS       *h2opus = (Mat_H2OPUS *)A->data;
1097036e4bdcSSatish Balay   PetscBool         isascii, vieweps;
1098f17f7ee2SSatish Balay   PetscMPIInt       size;
1099f17f7ee2SSatish Balay   PetscViewerFormat format;
1100f17f7ee2SSatish Balay 
1101f17f7ee2SSatish Balay   PetscFunctionBegin;
1102f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)view, PETSCVIEWERASCII, &isascii));
1103f17f7ee2SSatish Balay   PetscCall(PetscViewerGetFormat(view, &format));
1104f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1105f17f7ee2SSatish Balay   if (isascii) {
1106f17f7ee2SSatish Balay     if (format == PETSC_VIEWER_ASCII_MATLAB) {
1107f17f7ee2SSatish Balay       if (size == 1) {
1108f17f7ee2SSatish Balay         FILE *fp;
1109f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIGetPointer(view, &fp));
1110f17f7ee2SSatish Balay         dumpHMatrix(*h2opus->hmatrix, 6, fp);
1111f17f7ee2SSatish Balay       }
1112f17f7ee2SSatish Balay     } else {
1113f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  H-Matrix constructed from %s\n", h2opus->kernel ? "Kernel" : "Mat"));
1114f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  PointCloud dim %" PetscInt_FMT "\n", h2opus->ptcloud ? h2opus->ptcloud->getDimension() : 0));
1115f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  Admissibility parameters: leaf size %" PetscInt_FMT ", eta %g\n", h2opus->leafsize, (double)h2opus->eta));
1116f17f7ee2SSatish Balay       if (!h2opus->kernel) {
1117f17f7ee2SSatish Balay         PetscCall(PetscViewerASCIIPrintf(view, "  Sampling parameters: max_rank %" PetscInt_FMT ", samples %" PetscInt_FMT ", tolerance %g\n", h2opus->max_rank, h2opus->bs, (double)h2opus->rtol));
1118f17f7ee2SSatish Balay       } else {
11194cf0e950SBarry Smith         PetscCall(PetscViewerASCIIPrintf(view, "  Off-diagonal blocks approximation order %" PetscInt_FMT "\n", h2opus->basisord));
1120f17f7ee2SSatish Balay       }
1121f17f7ee2SSatish Balay       PetscCall(PetscViewerASCIIPrintf(view, "  Number of samples for norms %" PetscInt_FMT "\n", h2opus->norm_max_samples));
1122f17f7ee2SSatish Balay       if (size == 1) {
1123f17f7ee2SSatish Balay         double dense_mem_cpu = h2opus->hmatrix ? h2opus->hmatrix->getDenseMemoryUsage() : 0;
1124f17f7ee2SSatish Balay         double low_rank_cpu  = h2opus->hmatrix ? h2opus->hmatrix->getLowRankMemoryUsage() : 0;
1125f17f7ee2SSatish Balay   #if defined(PETSC_HAVE_CUDA)
1126f17f7ee2SSatish Balay         double dense_mem_gpu = h2opus->hmatrix_gpu ? h2opus->hmatrix_gpu->getDenseMemoryUsage() : 0;
1127f17f7ee2SSatish Balay         double low_rank_gpu  = h2opus->hmatrix_gpu ? h2opus->hmatrix_gpu->getLowRankMemoryUsage() : 0;
1128f17f7ee2SSatish Balay   #endif
1129f17f7ee2SSatish 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));
1130f17f7ee2SSatish Balay   #if defined(PETSC_HAVE_CUDA)
1131f17f7ee2SSatish 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));
1132f17f7ee2SSatish Balay   #endif
1133f17f7ee2SSatish Balay       } else {
1134f17f7ee2SSatish Balay   #if defined(PETSC_HAVE_CUDA)
1135f17f7ee2SSatish Balay         double      matrix_mem[4] = {0., 0., 0., 0.};
1136f17f7ee2SSatish Balay         PetscMPIInt rsize         = 4;
1137f17f7ee2SSatish Balay   #else
1138f17f7ee2SSatish Balay         double      matrix_mem[2] = {0., 0.};
1139f17f7ee2SSatish Balay         PetscMPIInt rsize         = 2;
1140f17f7ee2SSatish Balay   #endif
1141f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1142f17f7ee2SSatish Balay         matrix_mem[0] = h2opus->dist_hmatrix ? h2opus->dist_hmatrix->getLocalDenseMemoryUsage() : 0;
1143f17f7ee2SSatish Balay         matrix_mem[1] = h2opus->dist_hmatrix ? h2opus->dist_hmatrix->getLocalLowRankMemoryUsage() : 0;
1144f17f7ee2SSatish Balay     #if defined(PETSC_HAVE_CUDA)
1145f17f7ee2SSatish Balay         matrix_mem[2] = h2opus->dist_hmatrix_gpu ? h2opus->dist_hmatrix_gpu->getLocalDenseMemoryUsage() : 0;
1146f17f7ee2SSatish Balay         matrix_mem[3] = h2opus->dist_hmatrix_gpu ? h2opus->dist_hmatrix_gpu->getLocalLowRankMemoryUsage() : 0;
1147f17f7ee2SSatish Balay     #endif
1148f17f7ee2SSatish Balay   #endif
1149f17f7ee2SSatish Balay         PetscCall(MPIU_Allreduce(MPI_IN_PLACE, matrix_mem, rsize, MPI_DOUBLE_PRECISION, MPI_SUM, PetscObjectComm((PetscObject)A)));
1150f17f7ee2SSatish 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]));
1151f17f7ee2SSatish Balay   #if defined(PETSC_HAVE_CUDA)
1152f17f7ee2SSatish 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]));
1153f17f7ee2SSatish Balay   #endif
1154f17f7ee2SSatish Balay       }
1155f17f7ee2SSatish Balay     }
1156f17f7ee2SSatish Balay   }
1157036e4bdcSSatish Balay   vieweps = PETSC_FALSE;
1158036e4bdcSSatish Balay   PetscCall(PetscOptionsGetBool(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_vieweps", &vieweps, NULL));
1159036e4bdcSSatish Balay   if (vieweps) {
1160f17f7ee2SSatish Balay     char        filename[256];
1161f17f7ee2SSatish Balay     const char *name;
1162f17f7ee2SSatish Balay 
1163f17f7ee2SSatish Balay     PetscCall(PetscObjectGetName((PetscObject)A, &name));
1164f17f7ee2SSatish Balay     PetscCall(PetscSNPrintf(filename, sizeof(filename), "%s_structure.eps", name));
1165036e4bdcSSatish Balay     PetscCall(PetscOptionsGetString(((PetscObject)A)->options, ((PetscObject)A)->prefix, "-mat_h2opus_vieweps_filename", filename, sizeof(filename), NULL));
1166f17f7ee2SSatish Balay     outputEps(*h2opus->hmatrix, filename);
1167f17f7ee2SSatish Balay   }
11683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1169f17f7ee2SSatish Balay }
1170f17f7ee2SSatish Balay 
11718434afd1SBarry Smith static PetscErrorCode MatH2OpusSetCoords_H2OPUS(Mat A, PetscInt spacedim, const PetscReal coords[], PetscBool cdist, MatH2OpusKernelFn *kernel, void *kernelctx)
1172d71ae5a4SJacob Faibussowitsch {
1173f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus = (Mat_H2OPUS *)A->data;
1174f17f7ee2SSatish Balay   PetscReal  *gcoords;
1175f17f7ee2SSatish Balay   PetscInt    N;
1176f17f7ee2SSatish Balay   MPI_Comm    comm;
1177f17f7ee2SSatish Balay   PetscMPIInt size;
1178f17f7ee2SSatish Balay   PetscBool   cong;
1179f17f7ee2SSatish Balay 
1180f17f7ee2SSatish Balay   PetscFunctionBegin;
1181f17f7ee2SSatish Balay   PetscCall(PetscLayoutSetUp(A->rmap));
1182f17f7ee2SSatish Balay   PetscCall(PetscLayoutSetUp(A->cmap));
1183f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
1184f17f7ee2SSatish Balay   PetscCall(MatHasCongruentLayouts(A, &cong));
1185f17f7ee2SSatish Balay   PetscCheck(cong, comm, PETSC_ERR_SUP, "Only for square matrices with congruent layouts");
1186f17f7ee2SSatish Balay   N = A->rmap->N;
1187f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(comm, &size));
1188f17f7ee2SSatish Balay   if (spacedim > 0 && size > 1 && cdist) {
1189f17f7ee2SSatish Balay     PetscSF      sf;
1190f17f7ee2SSatish Balay     MPI_Datatype dtype;
1191f17f7ee2SSatish Balay 
1192f17f7ee2SSatish Balay     PetscCallMPI(MPI_Type_contiguous(spacedim, MPIU_REAL, &dtype));
1193f17f7ee2SSatish Balay     PetscCallMPI(MPI_Type_commit(&dtype));
1194f17f7ee2SSatish Balay 
1195f17f7ee2SSatish Balay     PetscCall(PetscSFCreate(comm, &sf));
1196f17f7ee2SSatish Balay     PetscCall(PetscSFSetGraphWithPattern(sf, A->rmap, PETSCSF_PATTERN_ALLGATHER));
1197f17f7ee2SSatish Balay     PetscCall(PetscMalloc1(spacedim * N, &gcoords));
1198f17f7ee2SSatish Balay     PetscCall(PetscSFBcastBegin(sf, dtype, coords, gcoords, MPI_REPLACE));
1199f17f7ee2SSatish Balay     PetscCall(PetscSFBcastEnd(sf, dtype, coords, gcoords, MPI_REPLACE));
1200f17f7ee2SSatish Balay     PetscCall(PetscSFDestroy(&sf));
1201f17f7ee2SSatish Balay     PetscCallMPI(MPI_Type_free(&dtype));
1202f17f7ee2SSatish Balay   } else gcoords = (PetscReal *)coords;
1203f17f7ee2SSatish Balay 
1204f17f7ee2SSatish Balay   delete h2opus->ptcloud;
1205f17f7ee2SSatish Balay   delete h2opus->kernel;
1206f17f7ee2SSatish Balay   h2opus->ptcloud = new PetscPointCloud<PetscReal>(spacedim, N, gcoords);
1207f17f7ee2SSatish Balay   if (kernel) h2opus->kernel = new PetscFunctionGenerator<PetscScalar>(kernel, spacedim, kernelctx);
1208f17f7ee2SSatish Balay   if (gcoords != coords) PetscCall(PetscFree(gcoords));
1209f17f7ee2SSatish Balay   A->preallocated = PETSC_TRUE;
12103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1211f17f7ee2SSatish Balay }
1212f17f7ee2SSatish Balay 
1213f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1214d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatBindToCPU_H2OPUS(Mat A, PetscBool flg)
1215d71ae5a4SJacob Faibussowitsch {
1216f17f7ee2SSatish Balay   PetscMPIInt size;
1217f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1218f17f7ee2SSatish Balay 
1219f17f7ee2SSatish Balay   PetscFunctionBegin;
1220f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1221f17f7ee2SSatish Balay   if (flg && A->offloadmask == PETSC_OFFLOAD_GPU) {
1222f17f7ee2SSatish Balay     if (size > 1) {
1223f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1224f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
1225f17f7ee2SSatish Balay       if (!a->dist_hmatrix) a->dist_hmatrix = new DistributedHMatrix(*a->dist_hmatrix_gpu);
1226f17f7ee2SSatish Balay       else *a->dist_hmatrix = *a->dist_hmatrix_gpu;
1227f17f7ee2SSatish Balay     #endif
1228f17f7ee2SSatish Balay     } else {
1229f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1230f17f7ee2SSatish Balay       if (!a->hmatrix) a->hmatrix = new HMatrix(*a->hmatrix_gpu);
1231f17f7ee2SSatish Balay       else *a->hmatrix = *a->hmatrix_gpu;
1232f17f7ee2SSatish Balay     }
1233f17f7ee2SSatish Balay     delete a->hmatrix_gpu;
1234f17f7ee2SSatish Balay     delete a->dist_hmatrix_gpu;
1235f17f7ee2SSatish Balay     a->hmatrix_gpu      = NULL;
1236f17f7ee2SSatish Balay     a->dist_hmatrix_gpu = NULL;
1237f17f7ee2SSatish Balay     A->offloadmask      = PETSC_OFFLOAD_CPU;
1238f17f7ee2SSatish Balay   } else if (!flg && A->offloadmask == PETSC_OFFLOAD_CPU) {
1239f17f7ee2SSatish Balay     if (size > 1) {
1240f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1241f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
1242f17f7ee2SSatish Balay       if (!a->dist_hmatrix_gpu) a->dist_hmatrix_gpu = new DistributedHMatrix_GPU(*a->dist_hmatrix);
1243f17f7ee2SSatish Balay       else *a->dist_hmatrix_gpu = *a->dist_hmatrix;
1244f17f7ee2SSatish Balay     #endif
1245f17f7ee2SSatish Balay     } else {
1246f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1247f17f7ee2SSatish Balay       if (!a->hmatrix_gpu) a->hmatrix_gpu = new HMatrix_GPU(*a->hmatrix);
1248f17f7ee2SSatish Balay       else *a->hmatrix_gpu = *a->hmatrix;
1249f17f7ee2SSatish Balay     }
1250f17f7ee2SSatish Balay     delete a->hmatrix;
1251f17f7ee2SSatish Balay     delete a->dist_hmatrix;
1252f17f7ee2SSatish Balay     a->hmatrix      = NULL;
1253f17f7ee2SSatish Balay     a->dist_hmatrix = NULL;
1254f17f7ee2SSatish Balay     A->offloadmask  = PETSC_OFFLOAD_GPU;
1255f17f7ee2SSatish Balay   }
1256f17f7ee2SSatish Balay   PetscCall(PetscFree(A->defaultvectype));
1257f17f7ee2SSatish Balay   if (!flg) {
1258f17f7ee2SSatish Balay     PetscCall(PetscStrallocpy(VECCUDA, &A->defaultvectype));
1259f17f7ee2SSatish Balay   } else {
1260f17f7ee2SSatish Balay     PetscCall(PetscStrallocpy(VECSTANDARD, &A->defaultvectype));
1261f17f7ee2SSatish Balay   }
1262f17f7ee2SSatish Balay   A->boundtocpu = flg;
12633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1264f17f7ee2SSatish Balay }
1265f17f7ee2SSatish Balay   #endif
1266f17f7ee2SSatish Balay 
1267f17f7ee2SSatish Balay /*MC
12681d27aa22SBarry Smith    MATH2OPUS = "h2opus" - A matrix type for hierarchical matrices using the H2Opus package {cite}`zampinibouakaramturkiyyahkniokeyes2022`.
1269f17f7ee2SSatish Balay 
12702ef1f0ffSBarry Smith    Options Database Key:
12712ef1f0ffSBarry Smith .  -mat_type h2opus - matrix type to "h2opus"
12722ef1f0ffSBarry Smith 
12732ef1f0ffSBarry Smith    Level: beginner
1274f17f7ee2SSatish Balay 
1275f17f7ee2SSatish Balay    Notes:
12761d27aa22SBarry Smith    H2Opus implements hierarchical matrices in the $H^2$ flavour. It supports CPU or NVIDIA GPUs.
127711a5261eSBarry Smith 
12782ef1f0ffSBarry Smith    For CPU only builds, use `./configure --download-h2opus --download-thrust` to install PETSc to use H2Opus.
12792ef1f0ffSBarry Smith    In order to run on NVIDIA GPUs, use `./configure --download-h2opus --download-magma --download-kblas`.
1280f17f7ee2SSatish Balay 
12811cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MATH2OPUS`, `MATHTOOL`, `MATDENSE`, `MatCreateH2OpusFromKernel()`, `MatCreateH2OpusFromMat()`
1282f17f7ee2SSatish Balay M*/
1283d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_H2OPUS(Mat A)
1284d71ae5a4SJacob Faibussowitsch {
1285f17f7ee2SSatish Balay   Mat_H2OPUS *a;
1286f17f7ee2SSatish Balay   PetscMPIInt size;
1287f17f7ee2SSatish Balay 
1288f17f7ee2SSatish Balay   PetscFunctionBegin;
1289f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1290f17f7ee2SSatish Balay   PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUDA));
1291f17f7ee2SSatish Balay   #endif
12924dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&a));
1293f17f7ee2SSatish Balay   A->data = (void *)a;
1294f17f7ee2SSatish Balay 
1295f17f7ee2SSatish Balay   a->eta              = 0.9;
1296f17f7ee2SSatish Balay   a->leafsize         = 32;
1297f17f7ee2SSatish Balay   a->basisord         = 4;
1298f17f7ee2SSatish Balay   a->max_rank         = 64;
1299f17f7ee2SSatish Balay   a->bs               = 32;
1300f17f7ee2SSatish Balay   a->rtol             = 1.e-4;
1301f17f7ee2SSatish Balay   a->s                = 1.0;
1302f17f7ee2SSatish Balay   a->norm_max_samples = 10;
1303036e4bdcSSatish Balay   a->resize           = PETSC_TRUE; /* reallocate after compression */
1304f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1305f17f7ee2SSatish Balay   h2opusCreateDistributedHandleComm(&a->handle, PetscObjectComm((PetscObject)A));
1306f17f7ee2SSatish Balay   #else
1307f17f7ee2SSatish Balay   h2opusCreateHandle(&a->handle);
1308f17f7ee2SSatish Balay   #endif
1309f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1310f17f7ee2SSatish Balay   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATH2OPUS));
1311f17f7ee2SSatish Balay   PetscCall(PetscMemzero(A->ops, sizeof(struct _MatOps)));
1312f17f7ee2SSatish Balay 
1313f17f7ee2SSatish Balay   A->ops->destroy          = MatDestroy_H2OPUS;
1314f17f7ee2SSatish Balay   A->ops->view             = MatView_H2OPUS;
1315f17f7ee2SSatish Balay   A->ops->assemblyend      = MatAssemblyEnd_H2OPUS;
1316f17f7ee2SSatish Balay   A->ops->mult             = MatMult_H2OPUS;
1317f17f7ee2SSatish Balay   A->ops->multtranspose    = MatMultTranspose_H2OPUS;
1318f17f7ee2SSatish Balay   A->ops->multadd          = MatMultAdd_H2OPUS;
1319f17f7ee2SSatish Balay   A->ops->multtransposeadd = MatMultTransposeAdd_H2OPUS;
1320f17f7ee2SSatish Balay   A->ops->scale            = MatScale_H2OPUS;
1321f17f7ee2SSatish Balay   A->ops->duplicate        = MatDuplicate_H2OPUS;
1322f17f7ee2SSatish Balay   A->ops->setfromoptions   = MatSetFromOptions_H2OPUS;
1323f17f7ee2SSatish Balay   A->ops->norm             = MatNorm_H2OPUS;
1324f17f7ee2SSatish Balay   A->ops->zeroentries      = MatZeroEntries_H2OPUS;
1325f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1326f17f7ee2SSatish Balay   A->ops->bindtocpu = MatBindToCPU_H2OPUS;
1327f17f7ee2SSatish Balay   #endif
1328f17f7ee2SSatish Balay 
1329f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdense_C", MatProductSetFromOptions_H2OPUS));
1330f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_seqdensecuda_C", MatProductSetFromOptions_H2OPUS));
1331f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidense_C", MatProductSetFromOptions_H2OPUS));
1332f17f7ee2SSatish Balay   PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatProductSetFromOptions_h2opus_mpidensecuda_C", MatProductSetFromOptions_H2OPUS));
1333f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1334f17f7ee2SSatish Balay   PetscCall(PetscFree(A->defaultvectype));
1335f17f7ee2SSatish Balay   PetscCall(PetscStrallocpy(VECCUDA, &A->defaultvectype));
1336f17f7ee2SSatish Balay   #endif
13373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1338f17f7ee2SSatish Balay }
1339f17f7ee2SSatish Balay 
1340*cc4c1da9SBarry Smith /*@
1341f17f7ee2SSatish Balay   MatH2OpusOrthogonalize - Orthogonalize the basis tree of a hierarchical matrix.
1342f17f7ee2SSatish Balay 
1343f17f7ee2SSatish Balay   Input Parameter:
1344f17f7ee2SSatish Balay . A - the matrix
1345f17f7ee2SSatish Balay 
1346f17f7ee2SSatish Balay   Level: intermediate
1347f17f7ee2SSatish Balay 
13481cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusCompress()`
1349f17f7ee2SSatish Balay @*/
1350d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusOrthogonalize(Mat A)
1351d71ae5a4SJacob Faibussowitsch {
1352f17f7ee2SSatish Balay   PetscBool   ish2opus;
1353f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1354f17f7ee2SSatish Balay   PetscMPIInt size;
1355f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
1356f17f7ee2SSatish Balay 
1357f17f7ee2SSatish Balay   PetscFunctionBegin;
1358f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1359f17f7ee2SSatish Balay   PetscValidType(A, 1);
1360f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
13613ba16761SJacob Faibussowitsch   if (!ish2opus) PetscFunctionReturn(PETSC_SUCCESS);
13623ba16761SJacob Faibussowitsch   if (a->orthogonal) PetscFunctionReturn(PETSC_SUCCESS);
1363f17f7ee2SSatish Balay   HLibProfile::clear();
1364f17f7ee2SSatish Balay   PetscCall(PetscLogEventBegin(MAT_H2Opus_Orthog, A, 0, 0, 0));
1365f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1366f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
1367f17f7ee2SSatish Balay   #endif
1368f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1369f17f7ee2SSatish Balay   if (size > 1) {
1370f17f7ee2SSatish Balay     if (boundtocpu) {
1371f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1372f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1373f17f7ee2SSatish Balay       distributed_horthog(*a->dist_hmatrix, a->handle);
1374f17f7ee2SSatish Balay   #endif
1375f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1376f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1377f17f7ee2SSatish Balay     } else {
1378f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1379f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1380f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
1381f17f7ee2SSatish Balay       distributed_horthog(*a->dist_hmatrix_gpu, a->handle);
1382f17f7ee2SSatish Balay     #endif
1383f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1384f17f7ee2SSatish Balay   #endif
1385f17f7ee2SSatish Balay     }
1386f17f7ee2SSatish Balay   } else {
1387f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1388f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle->handle;
1389f17f7ee2SSatish Balay   #else
1390f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle;
1391f17f7ee2SSatish Balay   #endif
1392f17f7ee2SSatish Balay     if (boundtocpu) {
1393f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1394f17f7ee2SSatish Balay       horthog(*a->hmatrix, handle);
1395f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1396f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1397f17f7ee2SSatish Balay     } else {
1398f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1399f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1400f17f7ee2SSatish Balay       horthog(*a->hmatrix_gpu, handle);
1401f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1402f17f7ee2SSatish Balay   #endif
1403f17f7ee2SSatish Balay     }
1404f17f7ee2SSatish Balay   }
1405f17f7ee2SSatish Balay   a->orthogonal = PETSC_TRUE;
1406f17f7ee2SSatish Balay   { /* log flops */
1407f17f7ee2SSatish Balay     double gops, time, perf, dev;
1408f17f7ee2SSatish Balay     HLibProfile::getHorthogPerf(gops, time, perf, dev);
1409f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1410f17f7ee2SSatish Balay     if (boundtocpu) {
1411f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
1412f17f7ee2SSatish Balay     } else {
1413f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
1414f17f7ee2SSatish Balay     }
1415f17f7ee2SSatish Balay   #else
1416f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
1417f17f7ee2SSatish Balay   #endif
1418f17f7ee2SSatish Balay   }
1419f17f7ee2SSatish Balay   PetscCall(PetscLogEventEnd(MAT_H2Opus_Orthog, A, 0, 0, 0));
14203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1421f17f7ee2SSatish Balay }
1422f17f7ee2SSatish Balay 
1423*cc4c1da9SBarry Smith /*@
1424f17f7ee2SSatish Balay   MatH2OpusCompress - Compress a hierarchical matrix.
1425f17f7ee2SSatish Balay 
1426f17f7ee2SSatish Balay   Input Parameters:
1427f17f7ee2SSatish Balay + A   - the matrix
1428f17f7ee2SSatish Balay - tol - the absolute truncation threshold
1429f17f7ee2SSatish Balay 
1430f17f7ee2SSatish Balay   Level: intermediate
1431f17f7ee2SSatish Balay 
14321cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusOrthogonalize()`
1433f17f7ee2SSatish Balay @*/
1434d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusCompress(Mat A, PetscReal tol)
1435d71ae5a4SJacob Faibussowitsch {
1436f17f7ee2SSatish Balay   PetscBool   ish2opus;
1437f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1438f17f7ee2SSatish Balay   PetscMPIInt size;
1439f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
1440f17f7ee2SSatish Balay 
1441f17f7ee2SSatish Balay   PetscFunctionBegin;
1442f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1443f17f7ee2SSatish Balay   PetscValidType(A, 1);
1444f17f7ee2SSatish Balay   PetscValidLogicalCollectiveReal(A, tol, 2);
1445f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
14463ba16761SJacob Faibussowitsch   if (!ish2opus || tol <= 0.0) PetscFunctionReturn(PETSC_SUCCESS);
1447f17f7ee2SSatish Balay   PetscCall(MatH2OpusOrthogonalize(A));
1448f17f7ee2SSatish Balay   HLibProfile::clear();
1449f17f7ee2SSatish Balay   PetscCall(PetscLogEventBegin(MAT_H2Opus_Compress, A, 0, 0, 0));
1450f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1451f17f7ee2SSatish Balay   boundtocpu = A->boundtocpu;
1452f17f7ee2SSatish Balay   #endif
1453f17f7ee2SSatish Balay   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1454f17f7ee2SSatish Balay   if (size > 1) {
1455f17f7ee2SSatish Balay     if (boundtocpu) {
1456f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1457f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1458f17f7ee2SSatish Balay       distributed_hcompress(*a->dist_hmatrix, tol, a->handle);
1459036e4bdcSSatish Balay       if (a->resize) {
1460036e4bdcSSatish Balay         DistributedHMatrix *dist_hmatrix = new DistributedHMatrix(*a->dist_hmatrix);
1461036e4bdcSSatish Balay         delete a->dist_hmatrix;
1462036e4bdcSSatish Balay         a->dist_hmatrix = dist_hmatrix;
1463036e4bdcSSatish Balay       }
1464f17f7ee2SSatish Balay   #endif
1465f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1466f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1467f17f7ee2SSatish Balay     } else {
1468f17f7ee2SSatish Balay       PetscCheck(a->dist_hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1469f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1470f17f7ee2SSatish Balay     #if defined(H2OPUS_USE_MPI)
1471f17f7ee2SSatish Balay       distributed_hcompress(*a->dist_hmatrix_gpu, tol, a->handle);
1472036e4bdcSSatish Balay 
1473036e4bdcSSatish Balay       if (a->resize) {
1474036e4bdcSSatish Balay         DistributedHMatrix_GPU *dist_hmatrix_gpu = new DistributedHMatrix_GPU(*a->dist_hmatrix_gpu);
1475036e4bdcSSatish Balay         delete a->dist_hmatrix_gpu;
1476036e4bdcSSatish Balay         a->dist_hmatrix_gpu = dist_hmatrix_gpu;
1477036e4bdcSSatish Balay       }
1478f17f7ee2SSatish Balay     #endif
1479f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1480f17f7ee2SSatish Balay   #endif
1481f17f7ee2SSatish Balay     }
1482f17f7ee2SSatish Balay   } else {
1483f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1484f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle->handle;
1485f17f7ee2SSatish Balay   #else
1486f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle;
1487f17f7ee2SSatish Balay   #endif
1488f17f7ee2SSatish Balay     if (boundtocpu) {
1489f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1490f17f7ee2SSatish Balay       hcompress(*a->hmatrix, tol, handle);
1491036e4bdcSSatish Balay 
1492036e4bdcSSatish Balay       if (a->resize) {
1493036e4bdcSSatish Balay         HMatrix *hmatrix = new HMatrix(*a->hmatrix);
1494036e4bdcSSatish Balay         delete a->hmatrix;
1495036e4bdcSSatish Balay         a->hmatrix = hmatrix;
1496036e4bdcSSatish Balay       }
1497f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1498f17f7ee2SSatish Balay       A->offloadmask = PETSC_OFFLOAD_CPU;
1499f17f7ee2SSatish Balay     } else {
1500f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1501f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeBegin());
1502f17f7ee2SSatish Balay       hcompress(*a->hmatrix_gpu, tol, handle);
1503f17f7ee2SSatish Balay       PetscCall(PetscLogGpuTimeEnd());
1504036e4bdcSSatish Balay 
1505036e4bdcSSatish Balay       if (a->resize) {
1506036e4bdcSSatish Balay         HMatrix_GPU *hmatrix_gpu = new HMatrix_GPU(*a->hmatrix_gpu);
1507036e4bdcSSatish Balay         delete a->hmatrix_gpu;
1508036e4bdcSSatish Balay         a->hmatrix_gpu = hmatrix_gpu;
1509036e4bdcSSatish Balay       }
1510f17f7ee2SSatish Balay   #endif
1511f17f7ee2SSatish Balay     }
1512f17f7ee2SSatish Balay   }
1513f17f7ee2SSatish Balay   { /* log flops */
1514f17f7ee2SSatish Balay     double gops, time, perf, dev;
1515f17f7ee2SSatish Balay     HLibProfile::getHcompressPerf(gops, time, perf, dev);
1516f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1517f17f7ee2SSatish Balay     if (boundtocpu) {
1518f17f7ee2SSatish Balay       PetscCall(PetscLogFlops(1e9 * gops));
1519f17f7ee2SSatish Balay     } else {
1520f17f7ee2SSatish Balay       PetscCall(PetscLogGpuFlops(1e9 * gops));
1521f17f7ee2SSatish Balay     }
1522f17f7ee2SSatish Balay   #else
1523f17f7ee2SSatish Balay     PetscCall(PetscLogFlops(1e9 * gops));
1524f17f7ee2SSatish Balay   #endif
1525f17f7ee2SSatish Balay   }
1526f17f7ee2SSatish Balay   PetscCall(PetscLogEventEnd(MAT_H2Opus_Compress, A, 0, 0, 0));
15273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1528f17f7ee2SSatish Balay }
1529f17f7ee2SSatish Balay 
1530*cc4c1da9SBarry Smith /*@
15312ef1f0ffSBarry Smith   MatH2OpusSetSamplingMat - Set a matrix to be sampled from matrix-vector products on another matrix to construct a hierarchical matrix.
1532f17f7ee2SSatish Balay 
1533f17f7ee2SSatish Balay   Input Parameters:
1534f17f7ee2SSatish Balay + A   - the hierarchical matrix
1535f17f7ee2SSatish Balay . B   - the matrix to be sampled
1536f17f7ee2SSatish Balay . bs  - maximum number of samples to be taken concurrently
1537f17f7ee2SSatish Balay - tol - relative tolerance for construction
1538f17f7ee2SSatish Balay 
15391d27aa22SBarry Smith   Level: intermediate
15401d27aa22SBarry Smith 
154111a5261eSBarry Smith   Notes:
154211a5261eSBarry Smith   You need to call `MatAssemblyBegin()` and `MatAssemblyEnd()` to update the hierarchical matrix.
1543f17f7ee2SSatish Balay 
15441cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusCompress()`, `MatH2OpusOrthogonalize()`
1545f17f7ee2SSatish Balay @*/
1546d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusSetSamplingMat(Mat A, Mat B, PetscInt bs, PetscReal tol)
1547d71ae5a4SJacob Faibussowitsch {
1548f17f7ee2SSatish Balay   PetscBool ish2opus;
1549f17f7ee2SSatish Balay 
1550f17f7ee2SSatish Balay   PetscFunctionBegin;
1551f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1552f17f7ee2SSatish Balay   PetscValidType(A, 1);
1553f17f7ee2SSatish Balay   if (B) PetscValidHeaderSpecific(B, MAT_CLASSID, 2);
1554f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(A, bs, 3);
155507e0ed13SBarry Smith   PetscValidLogicalCollectiveReal(A, tol, 4);
1556f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1557f17f7ee2SSatish Balay   if (ish2opus) {
1558f17f7ee2SSatish Balay     Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1559f17f7ee2SSatish Balay 
1560f17f7ee2SSatish Balay     if (!a->sampler) a->sampler = new PetscMatrixSampler();
1561f17f7ee2SSatish Balay     a->sampler->SetSamplingMat(B);
1562f17f7ee2SSatish Balay     if (bs > 0) a->bs = bs;
1563f17f7ee2SSatish Balay     if (tol > 0.) a->rtol = tol;
1564f17f7ee2SSatish Balay     delete a->kernel;
1565f17f7ee2SSatish Balay   }
15663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1567f17f7ee2SSatish Balay }
1568f17f7ee2SSatish Balay 
1569f17f7ee2SSatish Balay /*@C
157011a5261eSBarry Smith   MatCreateH2OpusFromKernel - Creates a `MATH2OPUS` from a user-supplied kernel.
1571f17f7ee2SSatish Balay 
1572f17f7ee2SSatish Balay   Input Parameters:
1573f17f7ee2SSatish Balay + comm      - MPI communicator
15742ef1f0ffSBarry Smith . m         - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
15752ef1f0ffSBarry Smith . n         - number of local columns (or `PETSC_DECIDE` to have calculated if `N` is given)
15762ef1f0ffSBarry Smith . M         - number of global rows (or `PETSC_DETERMINE` to have calculated if `m` is given)
15772ef1f0ffSBarry Smith . N         - number of global columns (or `PETSC_DETERMINE` to have calculated if `n` is given)
1578f17f7ee2SSatish Balay . spacedim  - dimension of the space coordinates
1579f17f7ee2SSatish Balay . coords    - coordinates of the points
1580f17f7ee2SSatish Balay . cdist     - whether or not coordinates are distributed
15812ef1f0ffSBarry Smith . kernel    - computational kernel (or `NULL`)
1582f17f7ee2SSatish Balay . kernelctx - kernel context
1583f17f7ee2SSatish Balay . eta       - admissibility condition tolerance
1584f17f7ee2SSatish Balay . leafsize  - leaf size in cluster tree
1585f17f7ee2SSatish Balay - basisord  - approximation order for Chebychev interpolation of low-rank blocks
1586f17f7ee2SSatish Balay 
1587f17f7ee2SSatish Balay   Output Parameter:
1588f17f7ee2SSatish Balay . nA - matrix
1589f17f7ee2SSatish Balay 
1590f17f7ee2SSatish Balay   Options Database Keys:
159111a5261eSBarry Smith + -mat_h2opus_leafsize <`PetscInt`>    - Leaf size of cluster tree
159211a5261eSBarry Smith . -mat_h2opus_eta <`PetscReal`>        - Admissibility condition tolerance
159311a5261eSBarry Smith . -mat_h2opus_order <`PetscInt`>       - Chebychev approximation order
159411a5261eSBarry Smith - -mat_h2opus_normsamples <`PetscInt`> - Maximum number of samples to be used when estimating norms
1595f17f7ee2SSatish Balay 
1596f17f7ee2SSatish Balay   Level: intermediate
1597f17f7ee2SSatish Balay 
15981cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`
1599f17f7ee2SSatish Balay @*/
16008434afd1SBarry Smith PetscErrorCode MatCreateH2OpusFromKernel(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt spacedim, const PetscReal coords[], PetscBool cdist, MatH2OpusKernelFn *kernel, void *kernelctx, PetscReal eta, PetscInt leafsize, PetscInt basisord, Mat *nA)
1601d71ae5a4SJacob Faibussowitsch {
1602f17f7ee2SSatish Balay   Mat         A;
1603f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus;
1604f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1605f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_FALSE;
1606f17f7ee2SSatish Balay   #else
1607f17f7ee2SSatish Balay   PetscBool iscpu = PETSC_TRUE;
1608f17f7ee2SSatish Balay   #endif
1609f17f7ee2SSatish Balay 
1610f17f7ee2SSatish Balay   PetscFunctionBegin;
1611036e4bdcSSatish Balay   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Different row and column local sizes are not supported");
1612f17f7ee2SSatish Balay   PetscCall(MatCreate(comm, &A));
1613f17f7ee2SSatish Balay   PetscCall(MatSetSizes(A, m, n, M, N));
1614036e4bdcSSatish Balay   PetscCheck(M == N, comm, PETSC_ERR_SUP, "Rectangular matrices are not supported");
1615f17f7ee2SSatish Balay   PetscCall(MatSetType(A, MATH2OPUS));
1616f17f7ee2SSatish Balay   PetscCall(MatBindToCPU(A, iscpu));
1617f17f7ee2SSatish Balay   PetscCall(MatH2OpusSetCoords_H2OPUS(A, spacedim, coords, cdist, kernel, kernelctx));
1618f17f7ee2SSatish Balay 
1619f17f7ee2SSatish Balay   h2opus = (Mat_H2OPUS *)A->data;
1620f17f7ee2SSatish Balay   if (eta > 0.) h2opus->eta = eta;
1621f17f7ee2SSatish Balay   if (leafsize > 0) h2opus->leafsize = leafsize;
1622f17f7ee2SSatish Balay   if (basisord > 0) h2opus->basisord = basisord;
1623f17f7ee2SSatish Balay 
1624f17f7ee2SSatish Balay   *nA = A;
16253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1626f17f7ee2SSatish Balay }
1627f17f7ee2SSatish Balay 
1628*cc4c1da9SBarry Smith /*@
162911a5261eSBarry Smith   MatCreateH2OpusFromMat - Creates a `MATH2OPUS` sampling from a user-supplied operator.
1630f17f7ee2SSatish Balay 
1631f17f7ee2SSatish Balay   Input Parameters:
1632f17f7ee2SSatish Balay + B        - the matrix to be sampled
1633f17f7ee2SSatish Balay . spacedim - dimension of the space coordinates
1634f17f7ee2SSatish Balay . coords   - coordinates of the points
1635f17f7ee2SSatish Balay . cdist    - whether or not coordinates are distributed
1636f17f7ee2SSatish Balay . eta      - admissibility condition tolerance
1637f17f7ee2SSatish Balay . leafsize - leaf size in cluster tree
1638f17f7ee2SSatish Balay . maxrank  - maximum rank allowed
1639f17f7ee2SSatish Balay . bs       - maximum number of samples to be taken concurrently
1640f17f7ee2SSatish Balay - rtol     - relative tolerance for construction
1641f17f7ee2SSatish Balay 
1642f17f7ee2SSatish Balay   Output Parameter:
1643f17f7ee2SSatish Balay . nA - matrix
1644f17f7ee2SSatish Balay 
1645f17f7ee2SSatish Balay   Options Database Keys:
164611a5261eSBarry Smith + -mat_h2opus_leafsize <`PetscInt`>      - Leaf size of cluster tree
164711a5261eSBarry Smith . -mat_h2opus_eta <`PetscReal`>          - Admissibility condition tolerance
164811a5261eSBarry Smith . -mat_h2opus_maxrank <`PetscInt`>       - Maximum rank when constructed from matvecs
164911a5261eSBarry Smith . -mat_h2opus_samples <`PetscInt`>       - Maximum number of samples to be taken concurrently when constructing from matvecs
165011a5261eSBarry Smith . -mat_h2opus_rtol <`PetscReal`>         - Relative tolerance for construction from sampling
165111a5261eSBarry Smith . -mat_h2opus_check <`PetscBool`>        - Check error when constructing from sampling during MatAssemblyEnd()
165211a5261eSBarry Smith . -mat_h2opus_hara_verbose <`PetscBool`> - Verbose output from hara construction
1653aaa8cc7dSPierre Jolivet - -mat_h2opus_normsamples <`PetscInt`>   - Maximum number of samples to be when estimating norms
1654f17f7ee2SSatish Balay 
16552ef1f0ffSBarry Smith   Level: intermediate
16562ef1f0ffSBarry Smith 
165711a5261eSBarry Smith   Note:
165811a5261eSBarry Smith   Not available in parallel
1659f17f7ee2SSatish Balay 
16601cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromKernel()`
1661f17f7ee2SSatish Balay @*/
1662d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateH2OpusFromMat(Mat B, PetscInt spacedim, const PetscReal coords[], PetscBool cdist, PetscReal eta, PetscInt leafsize, PetscInt maxrank, PetscInt bs, PetscReal rtol, Mat *nA)
1663d71ae5a4SJacob Faibussowitsch {
1664f17f7ee2SSatish Balay   Mat         A;
1665f17f7ee2SSatish Balay   Mat_H2OPUS *h2opus;
1666f17f7ee2SSatish Balay   MPI_Comm    comm;
1667f17f7ee2SSatish Balay   PetscBool   boundtocpu = PETSC_TRUE;
1668f17f7ee2SSatish Balay 
1669f17f7ee2SSatish Balay   PetscFunctionBegin;
1670f17f7ee2SSatish Balay   PetscValidHeaderSpecific(B, MAT_CLASSID, 1);
1671f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, spacedim, 2);
1672f17f7ee2SSatish Balay   PetscValidLogicalCollectiveBool(B, cdist, 4);
1673f17f7ee2SSatish Balay   PetscValidLogicalCollectiveReal(B, eta, 5);
1674f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, leafsize, 6);
1675f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, maxrank, 7);
1676f17f7ee2SSatish Balay   PetscValidLogicalCollectiveInt(B, bs, 8);
1677f17f7ee2SSatish Balay   PetscValidLogicalCollectiveReal(B, rtol, 9);
16784f572ea9SToby Isaac   PetscAssertPointer(nA, 10);
1679f17f7ee2SSatish Balay   PetscCall(PetscObjectGetComm((PetscObject)B, &comm));
1680036e4bdcSSatish Balay   PetscCheck(B->rmap->n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Different row and column local sizes are not supported");
1681036e4bdcSSatish Balay   PetscCheck(B->rmap->N == B->cmap->N, comm, PETSC_ERR_SUP, "Rectangular matrices are not supported");
1682f17f7ee2SSatish Balay   PetscCall(MatCreate(comm, &A));
1683f17f7ee2SSatish Balay   PetscCall(MatSetSizes(A, B->rmap->n, B->cmap->n, B->rmap->N, B->cmap->N));
1684f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1685f17f7ee2SSatish Balay   {
1686f17f7ee2SSatish Balay     VecType   vtype;
16872592a002SStefano Zampini     PetscBool isstd, iscuda, iskok;
1688f17f7ee2SSatish Balay 
1689f17f7ee2SSatish Balay     PetscCall(MatGetVecType(B, &vtype));
16902592a002SStefano Zampini     PetscCall(PetscStrcmpAny(vtype, &isstd, VECSTANDARD, VECSEQ, VECMPI, ""));
16912592a002SStefano Zampini     PetscCall(PetscStrcmpAny(vtype, &iscuda, VECCUDA, VECSEQCUDA, VECMPICUDA, ""));
16922592a002SStefano Zampini     PetscCall(PetscStrcmpAny(vtype, &iskok, VECKOKKOS, VECSEQKOKKOS, VECMPIKOKKOS, ""));
16932592a002SStefano Zampini     PetscCheck(isstd || iscuda || iskok, comm, PETSC_ERR_SUP, "Not for type %s", vtype);
1694f17f7ee2SSatish Balay     if (iscuda && !B->boundtocpu) boundtocpu = PETSC_FALSE;
16952592a002SStefano Zampini     if (iskok && PetscDefined(HAVE_MACRO_KOKKOS_ENABLE_CUDA)) boundtocpu = PETSC_FALSE;
1696f17f7ee2SSatish Balay   }
1697f17f7ee2SSatish Balay   #endif
1698f17f7ee2SSatish Balay   PetscCall(MatSetType(A, MATH2OPUS));
1699f17f7ee2SSatish Balay   PetscCall(MatBindToCPU(A, boundtocpu));
170048a46eb9SPierre Jolivet   if (spacedim) PetscCall(MatH2OpusSetCoords_H2OPUS(A, spacedim, coords, cdist, NULL, NULL));
1701f17f7ee2SSatish Balay   PetscCall(MatPropagateSymmetryOptions(B, A));
1702f17f7ee2SSatish Balay   /* PetscCheck(A->symmetric,comm,PETSC_ERR_SUP,"Unsymmetric sampling does not work"); */
1703f17f7ee2SSatish Balay 
1704f17f7ee2SSatish Balay   h2opus          = (Mat_H2OPUS *)A->data;
1705f17f7ee2SSatish Balay   h2opus->sampler = new PetscMatrixSampler(B);
1706f17f7ee2SSatish Balay   if (eta > 0.) h2opus->eta = eta;
1707f17f7ee2SSatish Balay   if (leafsize > 0) h2opus->leafsize = leafsize;
1708f17f7ee2SSatish Balay   if (maxrank > 0) h2opus->max_rank = maxrank;
1709f17f7ee2SSatish Balay   if (bs > 0) h2opus->bs = bs;
1710f17f7ee2SSatish Balay   if (rtol > 0.) h2opus->rtol = rtol;
1711f17f7ee2SSatish Balay   *nA             = A;
1712f17f7ee2SSatish Balay   A->preallocated = PETSC_TRUE;
17133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1714f17f7ee2SSatish Balay }
1715f17f7ee2SSatish Balay 
1716*cc4c1da9SBarry Smith /*@
1717f17f7ee2SSatish Balay   MatH2OpusGetIndexMap - Access reordering index set.
1718f17f7ee2SSatish Balay 
17192fe279fdSBarry Smith   Input Parameter:
1720f17f7ee2SSatish Balay . A - the matrix
1721f17f7ee2SSatish Balay 
1722f17f7ee2SSatish Balay   Output Parameter:
1723f17f7ee2SSatish Balay . indexmap - the index set for the reordering
1724f17f7ee2SSatish Balay 
1725f17f7ee2SSatish Balay   Level: intermediate
1726f17f7ee2SSatish Balay 
17271cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`
1728f17f7ee2SSatish Balay @*/
1729d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusGetIndexMap(Mat A, IS *indexmap)
1730d71ae5a4SJacob Faibussowitsch {
1731f17f7ee2SSatish Balay   PetscBool   ish2opus;
1732f17f7ee2SSatish Balay   Mat_H2OPUS *a = (Mat_H2OPUS *)A->data;
1733f17f7ee2SSatish Balay 
1734f17f7ee2SSatish Balay   PetscFunctionBegin;
1735f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1736f17f7ee2SSatish Balay   PetscValidType(A, 1);
17374f572ea9SToby Isaac   PetscAssertPointer(indexmap, 2);
1738f17f7ee2SSatish Balay   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1739f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1740f17f7ee2SSatish Balay   PetscCheck(ish2opus, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name);
1741f17f7ee2SSatish Balay   *indexmap = a->h2opus_indexmap;
17423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1743f17f7ee2SSatish Balay }
1744f17f7ee2SSatish Balay 
1745*cc4c1da9SBarry Smith /*@
1746f17f7ee2SSatish Balay   MatH2OpusMapVec - Maps a vector between PETSc and H2Opus ordering
1747f17f7ee2SSatish Balay 
1748f17f7ee2SSatish Balay   Input Parameters:
1749f17f7ee2SSatish Balay + A             - the matrix
1750f17f7ee2SSatish Balay . nativetopetsc - if true, maps from H2Opus ordering to PETSc ordering. If false, applies the reverse map
1751f17f7ee2SSatish Balay - in            - the vector to be mapped
1752f17f7ee2SSatish Balay 
1753f17f7ee2SSatish Balay   Output Parameter:
1754f17f7ee2SSatish Balay . out - the newly created mapped vector
1755f17f7ee2SSatish Balay 
1756f17f7ee2SSatish Balay   Level: intermediate
1757f17f7ee2SSatish Balay 
17581cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`
1759f17f7ee2SSatish Balay @*/
1760d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusMapVec(Mat A, PetscBool nativetopetsc, Vec in, Vec *out)
1761d71ae5a4SJacob Faibussowitsch {
1762f17f7ee2SSatish Balay   PetscBool    ish2opus;
1763f17f7ee2SSatish Balay   Mat_H2OPUS  *a = (Mat_H2OPUS *)A->data;
1764f17f7ee2SSatish Balay   PetscScalar *xin, *xout;
1765f17f7ee2SSatish Balay   PetscBool    nm;
1766f17f7ee2SSatish Balay 
1767f17f7ee2SSatish Balay   PetscFunctionBegin;
1768f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1769f17f7ee2SSatish Balay   PetscValidType(A, 1);
1770f17f7ee2SSatish Balay   PetscValidLogicalCollectiveBool(A, nativetopetsc, 2);
1771f17f7ee2SSatish Balay   PetscValidHeaderSpecific(in, VEC_CLASSID, 3);
17724f572ea9SToby Isaac   PetscAssertPointer(out, 4);
1773f17f7ee2SSatish Balay   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1774f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &ish2opus));
1775f17f7ee2SSatish Balay   PetscCheck(ish2opus, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not for type %s", ((PetscObject)A)->type_name);
1776f17f7ee2SSatish Balay   nm = a->nativemult;
1777f17f7ee2SSatish Balay   PetscCall(MatH2OpusSetNativeMult(A, (PetscBool)!nativetopetsc));
1778f17f7ee2SSatish Balay   PetscCall(MatCreateVecs(A, out, NULL));
1779f17f7ee2SSatish Balay   PetscCall(MatH2OpusSetNativeMult(A, nm));
1780f17f7ee2SSatish Balay   if (!a->sf) { /* same ordering */
1781f17f7ee2SSatish Balay     PetscCall(VecCopy(in, *out));
17823ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1783f17f7ee2SSatish Balay   }
1784f17f7ee2SSatish Balay   PetscCall(VecGetArrayRead(in, (const PetscScalar **)&xin));
1785f17f7ee2SSatish Balay   PetscCall(VecGetArrayWrite(*out, &xout));
1786f17f7ee2SSatish Balay   if (nativetopetsc) {
1787f17f7ee2SSatish Balay     PetscCall(PetscSFReduceBegin(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1788f17f7ee2SSatish Balay     PetscCall(PetscSFReduceEnd(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1789f17f7ee2SSatish Balay   } else {
1790f17f7ee2SSatish Balay     PetscCall(PetscSFBcastBegin(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1791f17f7ee2SSatish Balay     PetscCall(PetscSFBcastEnd(a->sf, MPIU_SCALAR, xin, xout, MPI_REPLACE));
1792f17f7ee2SSatish Balay   }
1793f17f7ee2SSatish Balay   PetscCall(VecRestoreArrayRead(in, (const PetscScalar **)&xin));
1794f17f7ee2SSatish Balay   PetscCall(VecRestoreArrayWrite(*out, &xout));
17953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1796f17f7ee2SSatish Balay }
1797f17f7ee2SSatish Balay 
1798*cc4c1da9SBarry Smith /*@
17991d27aa22SBarry Smith   MatH2OpusLowRankUpdate - Perform a low-rank update of the form $ A = A + s * U * V^T $
1800f17f7ee2SSatish Balay 
1801f17f7ee2SSatish Balay   Input Parameters:
180211a5261eSBarry Smith + A - the hierarchical `MATH2OPUS` matrix
1803f17f7ee2SSatish Balay . s - the scaling factor
1804f17f7ee2SSatish Balay . U - the dense low-rank update matrix
18052ef1f0ffSBarry Smith - V - (optional) the dense low-rank update matrix (if `NULL`, then `V` = `U` is assumed)
1806f17f7ee2SSatish Balay 
180711a5261eSBarry Smith   Note:
18081d27aa22SBarry Smith   The `U` and `V` matrices must be in `MATDENSE` dense format
1809f17f7ee2SSatish Balay 
1810f17f7ee2SSatish Balay   Level: intermediate
1811f17f7ee2SSatish Balay 
18121cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatCreate()`, `MATH2OPUS`, `MatCreateH2OpusFromMat()`, `MatCreateH2OpusFromKernel()`, `MatH2OpusCompress()`, `MatH2OpusOrthogonalize()`, `MATDENSE`
1813f17f7ee2SSatish Balay @*/
1814d71ae5a4SJacob Faibussowitsch PetscErrorCode MatH2OpusLowRankUpdate(Mat A, Mat U, Mat V, PetscScalar s)
1815d71ae5a4SJacob Faibussowitsch {
1816f17f7ee2SSatish Balay   PetscBool flg;
1817f17f7ee2SSatish Balay 
1818f17f7ee2SSatish Balay   PetscFunctionBegin;
1819f17f7ee2SSatish Balay   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
1820f17f7ee2SSatish Balay   PetscValidType(A, 1);
1821f17f7ee2SSatish Balay   PetscCheck(A->assembled, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
1822f17f7ee2SSatish Balay   PetscValidHeaderSpecific(U, MAT_CLASSID, 2);
1823f17f7ee2SSatish Balay   PetscCheckSameComm(A, 1, U, 2);
1824f17f7ee2SSatish Balay   if (V) {
1825f17f7ee2SSatish Balay     PetscValidHeaderSpecific(V, MAT_CLASSID, 3);
1826f17f7ee2SSatish Balay     PetscCheckSameComm(A, 1, V, 3);
1827f17f7ee2SSatish Balay   }
1828f17f7ee2SSatish Balay   PetscValidLogicalCollectiveScalar(A, s, 4);
1829f17f7ee2SSatish Balay 
1830f17f7ee2SSatish Balay   if (!V) V = U;
1831036e4bdcSSatish 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);
18323ba16761SJacob Faibussowitsch   if (!U->cmap->N) PetscFunctionReturn(PETSC_SUCCESS);
1833f17f7ee2SSatish Balay   PetscCall(PetscLayoutCompare(U->rmap, A->rmap, &flg));
1834f17f7ee2SSatish Balay   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "A and U must have the same row layout");
1835f17f7ee2SSatish Balay   PetscCall(PetscLayoutCompare(V->rmap, A->cmap, &flg));
1836f17f7ee2SSatish Balay   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONGSTATE, "A column layout must match V row column layout");
1837f17f7ee2SSatish Balay   PetscCall(PetscObjectTypeCompare((PetscObject)A, MATH2OPUS, &flg));
1838f17f7ee2SSatish Balay   if (flg) {
1839f17f7ee2SSatish Balay     Mat_H2OPUS        *a = (Mat_H2OPUS *)A->data;
1840f17f7ee2SSatish Balay     const PetscScalar *u, *v, *uu, *vv;
1841f17f7ee2SSatish Balay     PetscInt           ldu, ldv;
1842f17f7ee2SSatish Balay     PetscMPIInt        size;
1843f17f7ee2SSatish Balay   #if defined(H2OPUS_USE_MPI)
1844f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle->handle;
1845f17f7ee2SSatish Balay   #else
1846f17f7ee2SSatish Balay     h2opusHandle_t handle = a->handle;
1847f17f7ee2SSatish Balay   #endif
1848f17f7ee2SSatish Balay     PetscBool usesf = (PetscBool)(a->sf && !a->nativemult);
1849f17f7ee2SSatish Balay     PetscSF   usf, vsf;
1850f17f7ee2SSatish Balay 
1851f17f7ee2SSatish Balay     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A), &size));
1852036e4bdcSSatish Balay     PetscCheck(size <= 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Not yet implemented in parallel");
1853f17f7ee2SSatish Balay     PetscCall(PetscLogEventBegin(MAT_H2Opus_LR, A, 0, 0, 0));
1854f17f7ee2SSatish Balay     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)U, &flg, MATSEQDENSE, MATMPIDENSE, ""));
1855f17f7ee2SSatish Balay     PetscCheck(flg, PetscObjectComm((PetscObject)U), PETSC_ERR_SUP, "Not for U of type %s", ((PetscObject)U)->type_name);
1856f17f7ee2SSatish Balay     PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)V, &flg, MATSEQDENSE, MATMPIDENSE, ""));
1857f17f7ee2SSatish Balay     PetscCheck(flg, PetscObjectComm((PetscObject)V), PETSC_ERR_SUP, "Not for V of type %s", ((PetscObject)V)->type_name);
1858f17f7ee2SSatish Balay     PetscCall(MatDenseGetLDA(U, &ldu));
1859f17f7ee2SSatish Balay     PetscCall(MatDenseGetLDA(V, &ldv));
1860f17f7ee2SSatish Balay     PetscCall(MatBoundToCPU(A, &flg));
1861f17f7ee2SSatish Balay     if (usesf) {
1862f17f7ee2SSatish Balay       PetscInt n;
1863f17f7ee2SSatish Balay 
18640dd791a8SStefano Zampini       PetscCall(MatDenseGetH2OpusStridedSF(U, a->sf, &usf));
18650dd791a8SStefano Zampini       PetscCall(MatDenseGetH2OpusStridedSF(V, a->sf, &vsf));
1866f17f7ee2SSatish Balay       PetscCall(MatH2OpusResizeBuffers_Private(A, U->cmap->N, V->cmap->N));
1867f17f7ee2SSatish Balay       PetscCall(PetscSFGetGraph(a->sf, NULL, &n, NULL, NULL));
1868f17f7ee2SSatish Balay       ldu = n;
1869f17f7ee2SSatish Balay       ldv = n;
1870f17f7ee2SSatish Balay     }
1871f17f7ee2SSatish Balay     if (flg) {
1872f17f7ee2SSatish Balay       PetscCheck(a->hmatrix, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing CPU matrix");
1873f17f7ee2SSatish Balay       PetscCall(MatDenseGetArrayRead(U, &u));
1874f17f7ee2SSatish Balay       PetscCall(MatDenseGetArrayRead(V, &v));
1875f17f7ee2SSatish Balay       if (usesf) {
1876f17f7ee2SSatish Balay         vv = MatH2OpusGetThrustPointer(*a->yy);
1877f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1878f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1879f17f7ee2SSatish Balay         if (U != V) {
1880f17f7ee2SSatish Balay           uu = MatH2OpusGetThrustPointer(*a->xx);
1881f17f7ee2SSatish Balay           PetscCall(PetscSFBcastBegin(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1882f17f7ee2SSatish Balay           PetscCall(PetscSFBcastEnd(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1883f17f7ee2SSatish Balay         } else uu = vv;
18849371c9d4SSatish Balay       } else {
18859371c9d4SSatish Balay         uu = u;
18869371c9d4SSatish Balay         vv = v;
18879371c9d4SSatish Balay       }
1888f17f7ee2SSatish Balay       hlru_global(*a->hmatrix, uu, ldu, vv, ldv, U->cmap->N, s, handle);
1889f17f7ee2SSatish Balay       PetscCall(MatDenseRestoreArrayRead(U, &u));
1890f17f7ee2SSatish Balay       PetscCall(MatDenseRestoreArrayRead(V, &v));
1891f17f7ee2SSatish Balay     } else {
1892f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1893f17f7ee2SSatish Balay       PetscBool flgU, flgV;
1894f17f7ee2SSatish Balay 
1895f17f7ee2SSatish Balay       PetscCheck(a->hmatrix_gpu, PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "Missing GPU matrix");
1896f17f7ee2SSatish Balay       PetscCall(PetscObjectTypeCompareAny((PetscObject)U, &flgU, MATSEQDENSE, MATMPIDENSE, ""));
1897f17f7ee2SSatish Balay       if (flgU) PetscCall(MatConvert(U, MATDENSECUDA, MAT_INPLACE_MATRIX, &U));
1898f17f7ee2SSatish Balay       PetscCall(PetscObjectTypeCompareAny((PetscObject)V, &flgV, MATSEQDENSE, MATMPIDENSE, ""));
1899f17f7ee2SSatish Balay       if (flgV) PetscCall(MatConvert(V, MATDENSECUDA, MAT_INPLACE_MATRIX, &V));
1900f17f7ee2SSatish Balay       PetscCall(MatDenseCUDAGetArrayRead(U, &u));
1901f17f7ee2SSatish Balay       PetscCall(MatDenseCUDAGetArrayRead(V, &v));
1902f17f7ee2SSatish Balay       if (usesf) {
1903f17f7ee2SSatish Balay         vv = MatH2OpusGetThrustPointer(*a->yy_gpu);
1904f17f7ee2SSatish Balay         PetscCall(PetscSFBcastBegin(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1905f17f7ee2SSatish Balay         PetscCall(PetscSFBcastEnd(vsf, MPIU_SCALAR, v, (PetscScalar *)vv, MPI_REPLACE));
1906f17f7ee2SSatish Balay         if (U != V) {
1907f17f7ee2SSatish Balay           uu = MatH2OpusGetThrustPointer(*a->xx_gpu);
1908f17f7ee2SSatish Balay           PetscCall(PetscSFBcastBegin(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1909f17f7ee2SSatish Balay           PetscCall(PetscSFBcastEnd(usf, MPIU_SCALAR, u, (PetscScalar *)uu, MPI_REPLACE));
1910f17f7ee2SSatish Balay         } else uu = vv;
19119371c9d4SSatish Balay       } else {
19129371c9d4SSatish Balay         uu = u;
19139371c9d4SSatish Balay         vv = v;
19149371c9d4SSatish Balay       }
1915f17f7ee2SSatish Balay   #else
1916f17f7ee2SSatish Balay       SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_PLIB, "This should not happen");
1917f17f7ee2SSatish Balay   #endif
1918f17f7ee2SSatish Balay       hlru_global(*a->hmatrix_gpu, uu, ldu, vv, ldv, U->cmap->N, s, handle);
1919f17f7ee2SSatish Balay   #if defined(PETSC_H2OPUS_USE_GPU)
1920f17f7ee2SSatish Balay       PetscCall(MatDenseCUDARestoreArrayRead(U, &u));
1921f17f7ee2SSatish Balay       PetscCall(MatDenseCUDARestoreArrayRead(V, &v));
1922f17f7ee2SSatish Balay       if (flgU) PetscCall(MatConvert(U, MATDENSE, MAT_INPLACE_MATRIX, &U));
1923f17f7ee2SSatish Balay       if (flgV) PetscCall(MatConvert(V, MATDENSE, MAT_INPLACE_MATRIX, &V));
1924f17f7ee2SSatish Balay   #endif
1925f17f7ee2SSatish Balay     }
1926f17f7ee2SSatish Balay     PetscCall(PetscLogEventEnd(MAT_H2Opus_LR, A, 0, 0, 0));
1927f17f7ee2SSatish Balay     a->orthogonal = PETSC_FALSE;
1928f17f7ee2SSatish Balay   }
19293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1930f17f7ee2SSatish Balay }
1931f17f7ee2SSatish Balay #endif
1932