1d8f25ad8SToby Isaac const char help[] = "Tests PetscDTPTrimmedEvalJet()";
2d8f25ad8SToby Isaac
3d8f25ad8SToby Isaac #include <petscdt.h>
4d8f25ad8SToby Isaac #include <petscblaslapack.h>
5d8f25ad8SToby Isaac #include <petscmat.h>
6d8f25ad8SToby Isaac
constructTabulationAndMass(PetscInt dim,PetscInt deg,PetscInt form,PetscInt jetDegree,PetscInt npoints,const PetscReal * points,const PetscReal * weights,PetscInt * _Nb,PetscInt * _Nf,PetscInt * _Nk,PetscReal ** B,PetscScalar ** M)7d71ae5a4SJacob Faibussowitsch static PetscErrorCode constructTabulationAndMass(PetscInt dim, PetscInt deg, PetscInt form, PetscInt jetDegree, PetscInt npoints, const PetscReal *points, const PetscReal *weights, PetscInt *_Nb, PetscInt *_Nf, PetscInt *_Nk, PetscReal **B, PetscScalar **M)
8d71ae5a4SJacob Faibussowitsch {
9d8f25ad8SToby Isaac PetscInt Nf; // Number of form components
10d8f25ad8SToby Isaac PetscInt Nbpt; // number of trimmed polynomials
11d8f25ad8SToby Isaac PetscInt Nk; // jet size
12d8f25ad8SToby Isaac PetscReal *p_trimmed;
13d8f25ad8SToby Isaac
14d8f25ad8SToby Isaac PetscFunctionBegin;
159566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(form), &Nf));
169566063dSJacob Faibussowitsch PetscCall(PetscDTPTrimmedSize(dim, deg, form, &Nbpt));
179566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + jetDegree, dim, &Nk));
189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed));
199566063dSJacob Faibussowitsch PetscCall(PetscDTPTrimmedEvalJet(dim, npoints, points, deg, form, jetDegree, p_trimmed));
20d8f25ad8SToby Isaac
21d8f25ad8SToby Isaac // compute the direct mass matrix
22d8f25ad8SToby Isaac PetscScalar *M_trimmed;
239566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt * Nbpt, &M_trimmed));
24d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt; i++) {
25d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) {
26d8f25ad8SToby Isaac PetscReal v = 0.;
27d8f25ad8SToby Isaac
28d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) {
29d8f25ad8SToby Isaac const PetscReal *p_i = &p_trimmed[(i * Nf + f) * Nk * npoints];
30d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + f) * Nk * npoints];
31d8f25ad8SToby Isaac
32ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) v += p_i[pt] * p_j[pt] * weights[pt];
33d8f25ad8SToby Isaac }
34d8f25ad8SToby Isaac M_trimmed[i * Nbpt + j] += v;
35d8f25ad8SToby Isaac }
36d8f25ad8SToby Isaac }
37d8f25ad8SToby Isaac *_Nb = Nbpt;
38d8f25ad8SToby Isaac *_Nf = Nf;
39d8f25ad8SToby Isaac *_Nk = Nk;
40d8f25ad8SToby Isaac *B = p_trimmed;
41d8f25ad8SToby Isaac *M = M_trimmed;
423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
43d8f25ad8SToby Isaac }
44d8f25ad8SToby Isaac
test(PetscInt dim,PetscInt deg,PetscInt form,PetscInt jetDegree,PetscBool cond)45d71ae5a4SJacob Faibussowitsch static PetscErrorCode test(PetscInt dim, PetscInt deg, PetscInt form, PetscInt jetDegree, PetscBool cond)
46d71ae5a4SJacob Faibussowitsch {
47d8f25ad8SToby Isaac PetscQuadrature q;
48d8f25ad8SToby Isaac PetscInt npoints;
49d8f25ad8SToby Isaac const PetscReal *points;
50d8f25ad8SToby Isaac const PetscReal *weights;
51d8f25ad8SToby Isaac PetscInt Nf; // Number of form components
52d8f25ad8SToby Isaac PetscInt Nk; // jet size
53d8f25ad8SToby Isaac PetscInt Nbpt; // number of trimmed polynomials
54d8f25ad8SToby Isaac PetscReal *p_trimmed;
55d8f25ad8SToby Isaac PetscScalar *M_trimmed;
56d8f25ad8SToby Isaac PetscReal *p_scalar;
57d8f25ad8SToby Isaac PetscInt Nbp; // number of scalar polynomials
58d8f25ad8SToby Isaac PetscScalar *Mcopy;
59d8f25ad8SToby Isaac PetscScalar *M_moments;
60d8f25ad8SToby Isaac PetscReal frob_err = 0.;
61d8f25ad8SToby Isaac Mat mat_trimmed;
62d8f25ad8SToby Isaac Mat mat_moments_T;
63d8f25ad8SToby Isaac Mat AinvB;
64d8f25ad8SToby Isaac PetscInt Nbm1;
65d8f25ad8SToby Isaac Mat Mm1;
66d8f25ad8SToby Isaac PetscReal *p_trimmed_copy;
67d8f25ad8SToby Isaac PetscReal *M_moment_real;
68d8f25ad8SToby Isaac
69d8f25ad8SToby Isaac PetscFunctionBegin;
70d8f25ad8SToby Isaac // Construct an appropriate quadrature
719566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 2, -1., 1., &q));
729566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights));
73d8f25ad8SToby Isaac
749566063dSJacob Faibussowitsch PetscCall(constructTabulationAndMass(dim, deg, form, jetDegree, npoints, points, weights, &Nbpt, &Nf, &Nk, &p_trimmed, &M_trimmed));
75d8f25ad8SToby Isaac
769566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + deg, dim, &Nbp));
779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbp * Nk * npoints, &p_scalar));
789566063dSJacob Faibussowitsch PetscCall(PetscDTPKDEvalJet(dim, npoints, points, deg, jetDegree, p_scalar));
79d8f25ad8SToby Isaac
809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt * Nbpt, &Mcopy));
81d8f25ad8SToby Isaac // Print the condition numbers (useful for testing out different bases internally in PetscDTPTrimmedEvalJet())
82d8f25ad8SToby Isaac #if !defined(PETSC_USE_COMPLEX)
83d8f25ad8SToby Isaac if (cond) {
84d8f25ad8SToby Isaac PetscReal *S;
85d8f25ad8SToby Isaac PetscScalar *work;
866497c311SBarry Smith PetscBLASInt n, lwork, lierr;
87d8f25ad8SToby Isaac
886497c311SBarry Smith PetscCall(PetscBLASIntCast(Nbpt, &n));
896497c311SBarry Smith PetscCall(PetscBLASIntCast(5 * Nbpt, &lwork));
909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt, &S));
919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * Nbpt, &work));
929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt));
93d8f25ad8SToby Isaac
94792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("N", "N", &n, &n, Mcopy, &n, S, NULL, &n, NULL, &n, work, &lwork, &lierr));
95d8f25ad8SToby Isaac PetscReal cond = S[0] / S[Nbpt - 1];
9663a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", form %" PetscInt_FMT ": condition number %g\n", dim, deg, form, (double)cond));
979566063dSJacob Faibussowitsch PetscCall(PetscFree(work));
989566063dSJacob Faibussowitsch PetscCall(PetscFree(S));
99d8f25ad8SToby Isaac }
100d8f25ad8SToby Isaac #endif
101d8f25ad8SToby Isaac
102d8f25ad8SToby Isaac // compute the moments with the orthonormal polynomials
1039566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt * Nbp * Nf, &M_moments));
104d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbp; i++) {
105d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) {
106d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) {
107d8f25ad8SToby Isaac PetscReal v = 0.;
108d8f25ad8SToby Isaac const PetscReal *p_i = &p_scalar[i * Nk * npoints];
109d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + f) * Nk * npoints];
110d8f25ad8SToby Isaac
111ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) v += p_i[pt] * p_j[pt] * weights[pt];
112d8f25ad8SToby Isaac M_moments[(i * Nf + f) * Nbpt + j] += v;
113d8f25ad8SToby Isaac }
114d8f25ad8SToby Isaac }
115d8f25ad8SToby Isaac }
116d8f25ad8SToby Isaac
117d8f25ad8SToby Isaac // subtract M_moments^T * M_moments from M_trimmed: because the trimmed polynomials should be contained in
118d8f25ad8SToby Isaac // the full polynomials, the result should be zero
1199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt));
120d8f25ad8SToby Isaac {
1216497c311SBarry Smith PetscBLASInt m, n, k;
122d8f25ad8SToby Isaac PetscScalar mone = -1.;
123d8f25ad8SToby Isaac PetscScalar one = 1.;
124d8f25ad8SToby Isaac
1256497c311SBarry Smith PetscCall(PetscBLASIntCast(Nbpt, &m));
1266497c311SBarry Smith PetscCall(PetscBLASIntCast(Nbpt, &n));
1276497c311SBarry Smith PetscCall(PetscBLASIntCast(Nbp * Nf, &k));
128792fecdfSBarry Smith PetscCallBLAS("BLASgemm", BLASgemm_("N", "T", &m, &n, &k, &mone, M_moments, &m, M_moments, &m, &one, Mcopy, &m));
129d8f25ad8SToby Isaac }
130d8f25ad8SToby Isaac
131d8f25ad8SToby Isaac frob_err = 0.;
132d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt * Nbpt; i++) frob_err += PetscRealPart(Mcopy[i]) * PetscRealPart(Mcopy[i]);
133d8f25ad8SToby Isaac frob_err = PetscSqrtReal(frob_err);
134d8f25ad8SToby Isaac
1357a46b595SBarry Smith PetscCheck(frob_err <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", form %" PetscInt_FMT ": trimmed projection error %g", dim, deg, form, (double)frob_err);
136d8f25ad8SToby Isaac
137d8f25ad8SToby Isaac // P trimmed is also supposed to contain the polynomials of one degree less: construction M_moment[0:sub,:] * M_trimmed^{-1} * M_moments[0:sub,:]^T should be the identity matrix
1389566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt, M_trimmed, &mat_trimmed));
1399566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + deg - 1, dim, &Nbm1));
1409566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbm1 * Nf, M_moments, &mat_moments_T));
1419566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat_moments_T, MAT_DO_NOT_COPY_VALUES, &AinvB));
1429566063dSJacob Faibussowitsch PetscCall(MatLUFactor(mat_trimmed, NULL, NULL, NULL));
1439566063dSJacob Faibussowitsch PetscCall(MatMatSolve(mat_trimmed, mat_moments_T, AinvB));
144fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(mat_moments_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &Mm1));
1459566063dSJacob Faibussowitsch PetscCall(MatShift(Mm1, -1.));
1469566063dSJacob Faibussowitsch PetscCall(MatNorm(Mm1, NORM_FROBENIUS, &frob_err));
1477a46b595SBarry Smith PetscCheck(frob_err <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", form %" PetscInt_FMT ": trimmed reverse projection error %g", dim, deg, form, (double)frob_err);
1489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Mm1));
1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AinvB));
1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_moments_T));
151d8f25ad8SToby Isaac
152d8f25ad8SToby Isaac // The Koszul differential applied to P trimmed (Lambda k+1) should be contained in P trimmed (Lambda k)
153d8f25ad8SToby Isaac if (PetscAbsInt(form) < dim) {
154d8f25ad8SToby Isaac PetscInt Nf1, Nbpt1, Nk1;
155d8f25ad8SToby Isaac PetscReal *p_trimmed1;
156d8f25ad8SToby Isaac PetscScalar *M_trimmed1;
157d8f25ad8SToby Isaac PetscInt (*pattern)[3];
158d8f25ad8SToby Isaac PetscReal *p_koszul;
159d8f25ad8SToby Isaac PetscScalar *M_koszul;
160d8f25ad8SToby Isaac PetscScalar *M_k_moment;
161d8f25ad8SToby Isaac Mat mat_koszul;
162d8f25ad8SToby Isaac Mat mat_k_moment_T;
163d8f25ad8SToby Isaac Mat AinvB;
164d8f25ad8SToby Isaac Mat prod;
165d8f25ad8SToby Isaac
1669371c9d4SSatish Balay PetscCall(constructTabulationAndMass(dim, deg, form < 0 ? form - 1 : form + 1, 0, npoints, points, weights, &Nbpt1, &Nf1, &Nk1, &p_trimmed1, &M_trimmed1));
167d8f25ad8SToby Isaac
1689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nf1 * (PetscAbsInt(form) + 1), &pattern));
1699566063dSJacob Faibussowitsch PetscCall(PetscDTAltVInteriorPattern(dim, PetscAbsInt(form) + 1, pattern));
170d8f25ad8SToby Isaac
171d8f25ad8SToby Isaac // apply the Koszul operator
1729566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt1 * Nf * npoints, &p_koszul));
173d8f25ad8SToby Isaac for (PetscInt b = 0; b < Nbpt1; b++) {
174d8f25ad8SToby Isaac for (PetscInt a = 0; a < Nf1 * (PetscAbsInt(form) + 1); a++) {
175d8f25ad8SToby Isaac PetscInt i, j, k;
176d8f25ad8SToby Isaac PetscReal sign;
177d8f25ad8SToby Isaac PetscReal *p_i;
178d8f25ad8SToby Isaac const PetscReal *p_j;
179d8f25ad8SToby Isaac
180d8f25ad8SToby Isaac i = pattern[a][0];
181ad540459SPierre Jolivet if (form < 0) i = Nf - 1 - i;
182d8f25ad8SToby Isaac j = pattern[a][1];
183ad540459SPierre Jolivet if (form < 0) j = Nf1 - 1 - j;
184d8f25ad8SToby Isaac k = pattern[a][2] < 0 ? -(pattern[a][2] + 1) : pattern[a][2];
185d8f25ad8SToby Isaac sign = pattern[a][2] < 0 ? -1 : 1;
186ad540459SPierre Jolivet if (form < 0 && (i & 1) ^ (j & 1)) sign = -sign;
187d8f25ad8SToby Isaac
188d8f25ad8SToby Isaac p_i = &p_koszul[(b * Nf + i) * npoints];
189d8f25ad8SToby Isaac p_j = &p_trimmed1[(b * Nf1 + j) * npoints];
190ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) p_i[pt] += p_j[pt] * points[pt * dim + k] * sign;
191d8f25ad8SToby Isaac }
192d8f25ad8SToby Isaac }
193d8f25ad8SToby Isaac
194d8f25ad8SToby Isaac // mass matrix of the result
1959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt1 * Nbpt1, &M_koszul));
196d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt1; i++) {
197d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt1; j++) {
198d8f25ad8SToby Isaac PetscReal val = 0.;
199d8f25ad8SToby Isaac
200d8f25ad8SToby Isaac for (PetscInt v = 0; v < Nf; v++) {
201d8f25ad8SToby Isaac const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints];
202d8f25ad8SToby Isaac const PetscReal *p_j = &p_koszul[(j * Nf + v) * npoints];
203d8f25ad8SToby Isaac
204ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) val += p_i[pt] * p_j[pt] * weights[pt];
205d8f25ad8SToby Isaac }
206d8f25ad8SToby Isaac M_koszul[i * Nbpt1 + j] = val;
207d8f25ad8SToby Isaac }
208d8f25ad8SToby Isaac }
209d8f25ad8SToby Isaac
210d8f25ad8SToby Isaac // moment matrix between the result and P trimmed
2119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt * Nbpt1, &M_k_moment));
212d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt1; i++) {
213d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) {
214d8f25ad8SToby Isaac PetscReal val = 0.;
215d8f25ad8SToby Isaac
216d8f25ad8SToby Isaac for (PetscInt v = 0; v < Nf; v++) {
217d8f25ad8SToby Isaac const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints];
218d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + v) * Nk * npoints];
219d8f25ad8SToby Isaac
220ad540459SPierre Jolivet for (PetscInt pt = 0; pt < npoints; pt++) val += p_i[pt] * p_j[pt] * weights[pt];
221d8f25ad8SToby Isaac }
222d8f25ad8SToby Isaac M_k_moment[i * Nbpt + j] = val;
223d8f25ad8SToby Isaac }
224d8f25ad8SToby Isaac }
225d8f25ad8SToby Isaac
226d8f25ad8SToby Isaac // M_k_moment M_trimmed^{-1} M_k_moment^T == M_koszul
2279566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt1, Nbpt1, M_koszul, &mat_koszul));
2289566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt1, M_k_moment, &mat_k_moment_T));
2299566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat_k_moment_T, MAT_DO_NOT_COPY_VALUES, &AinvB));
2309566063dSJacob Faibussowitsch PetscCall(MatMatSolve(mat_trimmed, mat_k_moment_T, AinvB));
231fb842aefSJose E. Roman PetscCall(MatTransposeMatMult(mat_k_moment_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DETERMINE, &prod));
2329566063dSJacob Faibussowitsch PetscCall(MatAXPY(prod, -1., mat_koszul, SAME_NONZERO_PATTERN));
2339566063dSJacob Faibussowitsch PetscCall(MatNorm(prod, NORM_FROBENIUS, &frob_err));
234*966bd95aSPierre Jolivet PetscCheck(frob_err <= PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", forms (%" PetscInt_FMT ", %" PetscInt_FMT "): koszul projection error %g", dim, deg, form, form < 0 ? (form - 1) : (form + 1), (double)frob_err);
235d8f25ad8SToby Isaac
2369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&prod));
2379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AinvB));
2389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_k_moment_T));
2399566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_koszul));
2409566063dSJacob Faibussowitsch PetscCall(PetscFree(M_k_moment));
2419566063dSJacob Faibussowitsch PetscCall(PetscFree(M_koszul));
2429566063dSJacob Faibussowitsch PetscCall(PetscFree(p_koszul));
2439566063dSJacob Faibussowitsch PetscCall(PetscFree(pattern));
2449566063dSJacob Faibussowitsch PetscCall(PetscFree(p_trimmed1));
2459566063dSJacob Faibussowitsch PetscCall(PetscFree(M_trimmed1));
246d8f25ad8SToby Isaac }
247d8f25ad8SToby Isaac
248d8f25ad8SToby Isaac // M_moments has shape [Nbp][Nf][Nbpt]
249d8f25ad8SToby Isaac // p_scalar has shape [Nbp][Nk][npoints]
250d8f25ad8SToby Isaac // contracting on [Nbp] should be the same shape as
251d8f25ad8SToby Isaac // p_trimmed, which is [Nbpt][Nf][Nk][npoints]
2529566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed_copy));
2539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbp * Nf * Nbpt, &M_moment_real));
254ad540459SPierre Jolivet for (PetscInt i = 0; i < Nbp * Nf * Nbpt; i++) M_moment_real[i] = PetscRealPart(M_moments[i]);
255d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) {
2566497c311SBarry Smith PetscBLASInt m, n, k, lda, ldb, ldc;
257d8f25ad8SToby Isaac PetscReal alpha = 1.0;
258d8f25ad8SToby Isaac PetscReal beta = 1.0;
259d8f25ad8SToby Isaac
2606497c311SBarry Smith PetscCall(PetscBLASIntCast(Nk * npoints, &m));
2616497c311SBarry Smith PetscCall(PetscBLASIntCast(Nbpt, &n));
2626497c311SBarry Smith PetscCall(PetscBLASIntCast(Nbp, &k));
2636497c311SBarry Smith PetscCall(PetscBLASIntCast(Nk * npoints, &lda));
2646497c311SBarry Smith PetscCall(PetscBLASIntCast(Nf * Nbpt, &ldb));
2656497c311SBarry Smith PetscCall(PetscBLASIntCast(Nf * Nk * npoints, &ldc));
266792fecdfSBarry Smith PetscCallBLAS("BLASREALgemm", BLASREALgemm_("N", "T", &m, &n, &k, &alpha, p_scalar, &lda, &M_moment_real[f * Nbpt], &ldb, &beta, &p_trimmed_copy[f * Nk * npoints], &ldc));
267d8f25ad8SToby Isaac }
268d8f25ad8SToby Isaac frob_err = 0.;
269ad540459SPierre Jolivet for (PetscInt i = 0; i < Nbpt * Nf * Nk * npoints; i++) frob_err += (p_trimmed_copy[i] - p_trimmed[i]) * (p_trimmed_copy[i] - p_trimmed[i]);
270d8f25ad8SToby Isaac frob_err = PetscSqrtReal(frob_err);
271d8f25ad8SToby Isaac
2721baa6e33SBarry Smith PetscCheck(frob_err < 10 * PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", form %" PetscInt_FMT ": jet error %g", dim, deg, form, (double)frob_err);
273d8f25ad8SToby Isaac
2749566063dSJacob Faibussowitsch PetscCall(PetscFree(M_moment_real));
2759566063dSJacob Faibussowitsch PetscCall(PetscFree(p_trimmed_copy));
2769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_trimmed));
2779566063dSJacob Faibussowitsch PetscCall(PetscFree(Mcopy));
2789566063dSJacob Faibussowitsch PetscCall(PetscFree(M_moments));
2799566063dSJacob Faibussowitsch PetscCall(PetscFree(M_trimmed));
2809566063dSJacob Faibussowitsch PetscCall(PetscFree(p_trimmed));
2819566063dSJacob Faibussowitsch PetscCall(PetscFree(p_scalar));
2829566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q));
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
284d8f25ad8SToby Isaac }
285d8f25ad8SToby Isaac
main(int argc,char ** argv)286d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
287d71ae5a4SJacob Faibussowitsch {
288d8f25ad8SToby Isaac PetscInt max_dim = 3;
289d8f25ad8SToby Isaac PetscInt max_deg = 4;
290d8f25ad8SToby Isaac PetscInt k = 3;
291d8f25ad8SToby Isaac PetscBool cond = PETSC_FALSE;
292d8f25ad8SToby Isaac
293327415f7SBarry Smith PetscFunctionBeginUser;
2949566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
295d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for PetscDTPTrimmedEvalJet() tests", "none");
2969566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_dim", "Maximum dimension of the simplex", __FILE__, max_dim, &max_dim, NULL));
2979566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_degree", "Maximum degree of the trimmed polynomial space", __FILE__, max_deg, &max_deg, NULL));
2989566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_jet", "The number of derivatives to test", __FILE__, k, &k, NULL));
2999566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cond", "Compute the condition numbers of the mass matrices of the bases", __FILE__, cond, &cond, NULL));
300d0609cedSBarry Smith PetscOptionsEnd();
301d8f25ad8SToby Isaac for (PetscInt dim = 2; dim <= max_dim; dim++) {
302d8f25ad8SToby Isaac for (PetscInt deg = 1; deg <= max_deg; deg++) {
30348a46eb9SPierre Jolivet for (PetscInt form = -dim + 1; form <= dim; form++) PetscCall(test(dim, deg, form, PetscMax(1, k), cond));
304d8f25ad8SToby Isaac }
305d8f25ad8SToby Isaac }
3069566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
307b122ec5aSJacob Faibussowitsch return 0;
308d8f25ad8SToby Isaac }
309d8f25ad8SToby Isaac
310d8f25ad8SToby Isaac /*TEST
311d8f25ad8SToby Isaac
312d8f25ad8SToby Isaac test:
313d8f25ad8SToby Isaac requires: !single
314d8f25ad8SToby Isaac args:
3153886731fSPierre Jolivet output_file: output/empty.out
316d8f25ad8SToby Isaac
317d8f25ad8SToby Isaac TEST*/
318