xref: /petsc/src/dm/dt/tests/ex13.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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 
7d8f25ad8SToby Isaac static PetscErrorCode constructTabulationAndMass(PetscInt dim, PetscInt deg, PetscInt form, PetscInt jetDegree, PetscInt npoints,
8d8f25ad8SToby Isaac                                                  const PetscReal *points, const PetscReal *weights,
9d8f25ad8SToby Isaac                                                  PetscInt *_Nb, PetscInt *_Nf, PetscInt *_Nk,
10d8f25ad8SToby Isaac                                                  PetscReal **B, PetscScalar **M)
11d8f25ad8SToby Isaac {
12d8f25ad8SToby Isaac   PetscInt       Nf; // Number of form components
13d8f25ad8SToby Isaac   PetscInt       Nbpt; // number of trimmed polynomials
14d8f25ad8SToby Isaac   PetscInt       Nk; // jet size
15d8f25ad8SToby Isaac   PetscReal     *p_trimmed;
16d8f25ad8SToby Isaac 
17d8f25ad8SToby Isaac   PetscFunctionBegin;
18*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTBinomialInt(dim, PetscAbsInt(form), &Nf));
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTPTrimmedSize(dim, deg, form, &Nbpt));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTBinomialInt(dim + jetDegree, dim, &Nk));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTPTrimmedEvalJet(dim, npoints, points, deg, form, jetDegree, p_trimmed));
23d8f25ad8SToby Isaac 
24d8f25ad8SToby Isaac   // compute the direct mass matrix
25d8f25ad8SToby Isaac   PetscScalar *M_trimmed;
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(Nbpt * Nbpt, &M_trimmed));
27d8f25ad8SToby Isaac   for (PetscInt i = 0; i < Nbpt; i++) {
28d8f25ad8SToby Isaac     for (PetscInt j = 0; j < Nbpt; j++) {
29d8f25ad8SToby Isaac       PetscReal v = 0.;
30d8f25ad8SToby Isaac 
31d8f25ad8SToby Isaac       for (PetscInt f = 0; f < Nf; f++) {
32d8f25ad8SToby Isaac         const PetscReal *p_i = &p_trimmed[(i * Nf + f) * Nk * npoints];
33d8f25ad8SToby Isaac         const PetscReal *p_j = &p_trimmed[(j * Nf + f) * Nk * npoints];
34d8f25ad8SToby Isaac 
35d8f25ad8SToby Isaac         for (PetscInt pt = 0; pt < npoints; pt++) {
36d8f25ad8SToby Isaac           v += p_i[pt] * p_j[pt] * weights[pt];
37d8f25ad8SToby Isaac         }
38d8f25ad8SToby Isaac       }
39d8f25ad8SToby Isaac       M_trimmed[i * Nbpt + j] += v;
40d8f25ad8SToby Isaac     }
41d8f25ad8SToby Isaac   }
42d8f25ad8SToby Isaac   *_Nb = Nbpt;
43d8f25ad8SToby Isaac   *_Nf = Nf;
44d8f25ad8SToby Isaac   *_Nk = Nk;
45d8f25ad8SToby Isaac   *B = p_trimmed;
46d8f25ad8SToby Isaac   *M = M_trimmed;
47d8f25ad8SToby Isaac   PetscFunctionReturn(0);
48d8f25ad8SToby Isaac }
49d8f25ad8SToby Isaac 
50d8f25ad8SToby Isaac static PetscErrorCode test(PetscInt dim, PetscInt deg, PetscInt form, PetscInt jetDegree, PetscBool cond)
51d8f25ad8SToby Isaac {
52d8f25ad8SToby Isaac   PetscQuadrature  q;
53d8f25ad8SToby Isaac   PetscInt         npoints;
54d8f25ad8SToby Isaac   const PetscReal *points;
55d8f25ad8SToby Isaac   const PetscReal *weights;
56d8f25ad8SToby Isaac   PetscInt         Nf; // Number of form components
57d8f25ad8SToby Isaac   PetscInt         Nk; // jet size
58d8f25ad8SToby Isaac   PetscInt         Nbpt; // number of trimmed polynomials
59d8f25ad8SToby Isaac   PetscReal       *p_trimmed;
60d8f25ad8SToby Isaac   PetscScalar     *M_trimmed;
61d8f25ad8SToby Isaac   PetscReal       *p_scalar;
62d8f25ad8SToby Isaac   PetscInt         Nbp; // number of scalar polynomials
63d8f25ad8SToby Isaac   PetscScalar     *Mcopy;
64d8f25ad8SToby Isaac   PetscScalar     *M_moments;
65d8f25ad8SToby Isaac   PetscReal        frob_err = 0.;
66d8f25ad8SToby Isaac   Mat              mat_trimmed;
67d8f25ad8SToby Isaac   Mat              mat_moments_T;
68d8f25ad8SToby Isaac   Mat              AinvB;
69d8f25ad8SToby Isaac   PetscInt         Nbm1;
70d8f25ad8SToby Isaac   Mat              Mm1;
71d8f25ad8SToby Isaac   PetscReal       *p_trimmed_copy;
72d8f25ad8SToby Isaac   PetscReal       *M_moment_real;
73d8f25ad8SToby Isaac   PetscErrorCode   ierr;
74d8f25ad8SToby Isaac 
75d8f25ad8SToby Isaac   PetscFunctionBegin;
76d8f25ad8SToby Isaac   // Construct an appropriate quadrature
77*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTStroudConicalQuadrature(dim, 1, deg + 2, -1., 1., &q));
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights));
79d8f25ad8SToby Isaac 
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(constructTabulationAndMass(dim, deg, form, jetDegree, npoints, points, weights, &Nbpt, &Nf, &Nk, &p_trimmed, &M_trimmed));
81d8f25ad8SToby Isaac 
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTBinomialInt(dim + deg, dim, &Nbp));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(Nbp * Nk * npoints, &p_scalar));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTPKDEvalJet(dim, npoints, points, deg, jetDegree, p_scalar));
85d8f25ad8SToby Isaac 
86*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(Nbpt * Nbpt, &Mcopy));
87d8f25ad8SToby Isaac   // Print the condition numbers (useful for testing out different bases internally in PetscDTPTrimmedEvalJet())
88d8f25ad8SToby Isaac #if !defined(PETSC_USE_COMPLEX)
89d8f25ad8SToby Isaac   if (cond) {
90d8f25ad8SToby Isaac     PetscReal *S;
91d8f25ad8SToby Isaac     PetscScalar *work;
92d8f25ad8SToby Isaac     PetscBLASInt n = Nbpt;
93d8f25ad8SToby Isaac     PetscBLASInt lwork = 5 * Nbpt;
94d8f25ad8SToby Isaac     PetscBLASInt lierr;
95d8f25ad8SToby Isaac 
96*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(Nbpt, &S));
97*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(5*Nbpt, &work));
98*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt));
99d8f25ad8SToby Isaac 
100d8f25ad8SToby Isaac     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&n,&n,Mcopy,&n,S,NULL,&n,NULL,&n,work,&lwork,&lierr));
101d8f25ad8SToby Isaac     PetscReal cond = S[0] / S[Nbpt - 1];
102d8f25ad8SToby Isaac     ierr = PetscPrintf(PETSC_COMM_WORLD, "dimension %D, degree %D, form %D: condition number %g\n", dim, deg, form, (double) cond);
103*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(work));
104*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(S));
105d8f25ad8SToby Isaac   }
106d8f25ad8SToby Isaac #endif
107d8f25ad8SToby Isaac 
108d8f25ad8SToby Isaac   // compute the moments with the orthonormal polynomials
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(Nbpt * Nbp * Nf, &M_moments));
110d8f25ad8SToby Isaac   for (PetscInt i = 0; i < Nbp; i++) {
111d8f25ad8SToby Isaac     for (PetscInt j = 0; j < Nbpt; j++) {
112d8f25ad8SToby Isaac       for (PetscInt f = 0; f < Nf; f++) {
113d8f25ad8SToby Isaac         PetscReal        v = 0.;
114d8f25ad8SToby Isaac         const PetscReal *p_i = &p_scalar[i * Nk * npoints];
115d8f25ad8SToby Isaac         const PetscReal *p_j = &p_trimmed[(j * Nf + f) * Nk * npoints];
116d8f25ad8SToby Isaac 
117d8f25ad8SToby Isaac         for (PetscInt pt = 0; pt < npoints; pt++) {
118d8f25ad8SToby Isaac           v += p_i[pt] * p_j[pt] * weights[pt];
119d8f25ad8SToby Isaac         }
120d8f25ad8SToby Isaac         M_moments[(i * Nf + f) * Nbpt + j] += v;
121d8f25ad8SToby Isaac       }
122d8f25ad8SToby Isaac     }
123d8f25ad8SToby Isaac   }
124d8f25ad8SToby Isaac 
125d8f25ad8SToby Isaac   // subtract M_moments^T * M_moments from M_trimmed: because the trimmed polynomials should be contained in
126d8f25ad8SToby Isaac   // the full polynomials, the result should be zero
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt));
128d8f25ad8SToby Isaac   {
129d8f25ad8SToby Isaac     PetscBLASInt m = Nbpt;
130d8f25ad8SToby Isaac     PetscBLASInt n = Nbpt;
131d8f25ad8SToby Isaac     PetscBLASInt k = Nbp * Nf;
132d8f25ad8SToby Isaac     PetscScalar mone = -1.;
133d8f25ad8SToby Isaac     PetscScalar one = 1.;
134d8f25ad8SToby Isaac 
135d8f25ad8SToby Isaac     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&mone,M_moments,&m,M_moments,&m,&one,Mcopy,&m));
136d8f25ad8SToby Isaac   }
137d8f25ad8SToby Isaac 
138d8f25ad8SToby Isaac   frob_err = 0.;
139d8f25ad8SToby Isaac   for (PetscInt i = 0; i < Nbpt * Nbpt; i++) frob_err += PetscRealPart(Mcopy[i]) * PetscRealPart(Mcopy[i]);
140d8f25ad8SToby Isaac   frob_err = PetscSqrtReal(frob_err);
141d8f25ad8SToby Isaac 
142d8f25ad8SToby Isaac   if (frob_err > PETSC_SMALL) {
14398921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %D, degree %D, form %D: trimmed projection error %g", dim, deg, form, (double) frob_err);
144d8f25ad8SToby Isaac   }
145d8f25ad8SToby Isaac 
146d8f25ad8SToby 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
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt, M_trimmed, &mat_trimmed));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTBinomialInt(dim + deg - 1, dim, &Nbm1));
149*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbm1 * Nf, M_moments, &mat_moments_T));
150*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(mat_moments_T, MAT_DO_NOT_COPY_VALUES, &AinvB));
151*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLUFactor(mat_trimmed, NULL, NULL, NULL));
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMatSolve(mat_trimmed, mat_moments_T, AinvB));
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatTransposeMatMult(mat_moments_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Mm1));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(Mm1, -1.));
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(Mm1, NORM_FROBENIUS, &frob_err));
156d8f25ad8SToby Isaac   if (frob_err > PETSC_SMALL) {
15798921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %D, degree %D, form %D: trimmed reverse projection error %g", dim, deg, form, (double) frob_err);
158d8f25ad8SToby Isaac   }
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Mm1));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&AinvB));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat_moments_T));
162d8f25ad8SToby Isaac 
163d8f25ad8SToby Isaac   // The Koszul differential applied to P trimmed (Lambda k+1) should be contained in P trimmed (Lambda k)
164d8f25ad8SToby Isaac   if (PetscAbsInt(form) < dim) {
165d8f25ad8SToby Isaac     PetscInt     Nf1, Nbpt1, Nk1;
166d8f25ad8SToby Isaac     PetscReal   *p_trimmed1;
167d8f25ad8SToby Isaac     PetscScalar *M_trimmed1;
168d8f25ad8SToby Isaac     PetscInt   (*pattern)[3];
169d8f25ad8SToby Isaac     PetscReal   *p_koszul;
170d8f25ad8SToby Isaac     PetscScalar *M_koszul;
171d8f25ad8SToby Isaac     PetscScalar *M_k_moment;
172d8f25ad8SToby Isaac     Mat          mat_koszul;
173d8f25ad8SToby Isaac     Mat          mat_k_moment_T;
174d8f25ad8SToby Isaac     Mat          AinvB;
175d8f25ad8SToby Isaac     Mat          prod;
176d8f25ad8SToby Isaac 
177d8f25ad8SToby Isaac     ierr = constructTabulationAndMass(dim, deg, form < 0 ? form - 1 : form + 1, 0, npoints, points, weights, &Nbpt1, &Nf1, &Nk1,
178d8f25ad8SToby Isaac                                       &p_trimmed1, &M_trimmed1);CHKERRQ(ierr);
179d8f25ad8SToby Isaac 
180*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(Nf1 * (PetscAbsInt(form) + 1), &pattern));
181*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTAltVInteriorPattern(dim, PetscAbsInt(form) + 1, pattern));
182d8f25ad8SToby Isaac 
183d8f25ad8SToby Isaac     // apply the Koszul operator
184*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(Nbpt1 * Nf * npoints, &p_koszul));
185d8f25ad8SToby Isaac     for (PetscInt b = 0; b < Nbpt1; b++) {
186d8f25ad8SToby Isaac       for (PetscInt a = 0; a < Nf1 * (PetscAbsInt(form) + 1); a++) {
187d8f25ad8SToby Isaac         PetscInt         i,j,k;
188d8f25ad8SToby Isaac         PetscReal        sign;
189d8f25ad8SToby Isaac         PetscReal       *p_i;
190d8f25ad8SToby Isaac         const PetscReal *p_j;
191d8f25ad8SToby Isaac 
192d8f25ad8SToby Isaac         i = pattern[a][0];
193d8f25ad8SToby Isaac         if (form < 0) {
194d8f25ad8SToby Isaac           i = Nf-1-i;
195d8f25ad8SToby Isaac         }
196d8f25ad8SToby Isaac         j = pattern[a][1];
197d8f25ad8SToby Isaac         if (form < 0) {
198d8f25ad8SToby Isaac           j = Nf1-1-j;
199d8f25ad8SToby Isaac         }
200d8f25ad8SToby Isaac         k = pattern[a][2] < 0 ? -(pattern[a][2] + 1) : pattern[a][2];
201d8f25ad8SToby Isaac         sign = pattern[a][2] < 0 ? -1 : 1;
202d8f25ad8SToby Isaac         if (form < 0 && (i & 1) ^ (j & 1)) {
203d8f25ad8SToby Isaac           sign = -sign;
204d8f25ad8SToby Isaac         }
205d8f25ad8SToby Isaac 
206d8f25ad8SToby Isaac         p_i = &p_koszul[(b * Nf + i) * npoints];
207d8f25ad8SToby Isaac         p_j = &p_trimmed1[(b * Nf1 + j) * npoints];
208d8f25ad8SToby Isaac         for (PetscInt pt = 0; pt < npoints; pt++) {
209d8f25ad8SToby Isaac           p_i[pt] += p_j[pt] * points[pt * dim + k] * sign;
210d8f25ad8SToby Isaac         }
211d8f25ad8SToby Isaac       }
212d8f25ad8SToby Isaac     }
213d8f25ad8SToby Isaac 
214d8f25ad8SToby Isaac     // mass matrix of the result
215*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(Nbpt1 * Nbpt1, &M_koszul));
216d8f25ad8SToby Isaac     for (PetscInt i = 0; i < Nbpt1; i++) {
217d8f25ad8SToby Isaac       for (PetscInt j = 0; j < Nbpt1; j++) {
218d8f25ad8SToby Isaac         PetscReal val = 0.;
219d8f25ad8SToby Isaac 
220d8f25ad8SToby Isaac         for (PetscInt v = 0; v < Nf; v++) {
221d8f25ad8SToby Isaac           const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints];
222d8f25ad8SToby Isaac           const PetscReal *p_j = &p_koszul[(j * Nf + v) * npoints];
223d8f25ad8SToby Isaac 
224d8f25ad8SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
225d8f25ad8SToby Isaac             val += p_i[pt] * p_j[pt] * weights[pt];
226d8f25ad8SToby Isaac           }
227d8f25ad8SToby Isaac         }
228d8f25ad8SToby Isaac         M_koszul[i * Nbpt1 + j] = val;
229d8f25ad8SToby Isaac       }
230d8f25ad8SToby Isaac     }
231d8f25ad8SToby Isaac 
232d8f25ad8SToby Isaac     // moment matrix between the result and P trimmed
233*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(Nbpt * Nbpt1, &M_k_moment));
234d8f25ad8SToby Isaac     for (PetscInt i = 0; i < Nbpt1; i++) {
235d8f25ad8SToby Isaac       for (PetscInt j = 0; j < Nbpt; j++) {
236d8f25ad8SToby Isaac         PetscReal val = 0.;
237d8f25ad8SToby Isaac 
238d8f25ad8SToby Isaac         for (PetscInt v = 0; v < Nf; v++) {
239d8f25ad8SToby Isaac           const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints];
240d8f25ad8SToby Isaac           const PetscReal *p_j = &p_trimmed[(j * Nf + v) * Nk * npoints];
241d8f25ad8SToby Isaac 
242d8f25ad8SToby Isaac           for (PetscInt pt = 0; pt < npoints; pt++) {
243d8f25ad8SToby Isaac             val += p_i[pt] * p_j[pt] * weights[pt];
244d8f25ad8SToby Isaac           }
245d8f25ad8SToby Isaac         }
246d8f25ad8SToby Isaac         M_k_moment[i * Nbpt + j] = val;
247d8f25ad8SToby Isaac       }
248d8f25ad8SToby Isaac     }
249d8f25ad8SToby Isaac 
250d8f25ad8SToby Isaac     // M_k_moment M_trimmed^{-1} M_k_moment^T == M_koszul
251*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt1, Nbpt1, M_koszul, &mat_koszul));
252*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt1, M_k_moment, &mat_k_moment_T));
253*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDuplicate(mat_k_moment_T, MAT_DO_NOT_COPY_VALUES, &AinvB));
254*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMatSolve(mat_trimmed, mat_k_moment_T, AinvB));
255*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatTransposeMatMult(mat_k_moment_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &prod));
256*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAXPY(prod, -1., mat_koszul, SAME_NONZERO_PATTERN));
257*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatNorm(prod, NORM_FROBENIUS, &frob_err));
258d8f25ad8SToby Isaac     if (frob_err > PETSC_SMALL) {
25998921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %D, degree %D, forms (%D, %D): koszul projection error %g", dim, deg, form, form < 0 ? (form-1):(form+1), (double) frob_err);
260d8f25ad8SToby Isaac     }
261d8f25ad8SToby Isaac 
262*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&prod));
263*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&AinvB));
264*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&mat_k_moment_T));
265*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&mat_koszul));
266*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(M_k_moment));
267*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(M_koszul));
268*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(p_koszul));
269*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(pattern));
270*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(p_trimmed1));
271*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(M_trimmed1));
272d8f25ad8SToby Isaac   }
273d8f25ad8SToby Isaac 
274d8f25ad8SToby Isaac   // M_moments has shape [Nbp][Nf][Nbpt]
275d8f25ad8SToby Isaac   // p_scalar has shape [Nbp][Nk][npoints]
276d8f25ad8SToby Isaac   // contracting on [Nbp] should be the same shape as
277d8f25ad8SToby Isaac   // p_trimmed, which is [Nbpt][Nf][Nk][npoints]
278*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed_copy));
279*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(Nbp * Nf * Nbpt, &M_moment_real));
280d8f25ad8SToby Isaac   for (PetscInt i = 0; i < Nbp * Nf * Nbpt; i++) {
281d8f25ad8SToby Isaac     M_moment_real[i] = PetscRealPart(M_moments[i]);
282d8f25ad8SToby Isaac   }
283d8f25ad8SToby Isaac   for (PetscInt f = 0; f < Nf; f++) {
284d8f25ad8SToby Isaac     PetscBLASInt m = Nk * npoints;
285d8f25ad8SToby Isaac     PetscBLASInt n = Nbpt;
286d8f25ad8SToby Isaac     PetscBLASInt k = Nbp;
287d8f25ad8SToby Isaac     PetscBLASInt lda = Nk * npoints;
288d8f25ad8SToby Isaac     PetscBLASInt ldb = Nf * Nbpt;
289d8f25ad8SToby Isaac     PetscBLASInt ldc = Nf * Nk * npoints;
290d8f25ad8SToby Isaac     PetscReal    alpha = 1.0;
291d8f25ad8SToby Isaac     PetscReal    beta = 1.0;
292d8f25ad8SToby Isaac 
293d8f25ad8SToby Isaac     PetscStackCallBLAS("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));
294d8f25ad8SToby Isaac   }
295d8f25ad8SToby Isaac   frob_err = 0.;
296d8f25ad8SToby Isaac   for (PetscInt i = 0; i < Nbpt * Nf * Nk * npoints; i++) {
297d8f25ad8SToby Isaac     frob_err += (p_trimmed_copy[i] - p_trimmed[i]) * (p_trimmed_copy[i] - p_trimmed[i]);
298d8f25ad8SToby Isaac   }
299d8f25ad8SToby Isaac   frob_err = PetscSqrtReal(frob_err);
300d8f25ad8SToby Isaac 
301d8f25ad8SToby Isaac   if (frob_err > PETSC_SMALL) {
30298921bdaSJacob Faibussowitsch     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %D, degree %D, form %D: jet error %g", dim, deg, form, (double) frob_err);
303d8f25ad8SToby Isaac   }
304d8f25ad8SToby Isaac 
305*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(M_moment_real));
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(p_trimmed_copy));
307*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&mat_trimmed));
308*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(Mcopy));
309*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(M_moments));
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(M_trimmed));
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(p_trimmed));
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(p_scalar));
313*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscQuadratureDestroy(&q));
314d8f25ad8SToby Isaac   PetscFunctionReturn(0);
315d8f25ad8SToby Isaac }
316d8f25ad8SToby Isaac 
317d8f25ad8SToby Isaac int main(int argc, char **argv)
318d8f25ad8SToby Isaac {
319d8f25ad8SToby Isaac   PetscInt  max_dim = 3;
320d8f25ad8SToby Isaac   PetscInt  max_deg = 4;
321d8f25ad8SToby Isaac   PetscInt  k = 3;
322d8f25ad8SToby Isaac   PetscBool cond = PETSC_FALSE;
323d8f25ad8SToby Isaac 
324d8f25ad8SToby Isaac   PetscErrorCode ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
325d8f25ad8SToby Isaac   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for PetscDTPTrimmedEvalJet() tests","none");CHKERRQ(ierr);
326*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-max_dim", "Maximum dimension of the simplex",__FILE__,max_dim,&max_dim,NULL));
327*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-max_degree", "Maximum degree of the trimmed polynomial space",__FILE__,max_deg,&max_deg,NULL));
328*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsInt("-max_jet", "The number of derivatives to test",__FILE__,k,&k,NULL));
329*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-cond", "Compute the condition numbers of the mass matrices of the bases",__FILE__,cond,&cond,NULL));
330d8f25ad8SToby Isaac   ierr = PetscOptionsEnd();CHKERRQ(ierr);
331d8f25ad8SToby Isaac   for (PetscInt dim = 2; dim <= max_dim; dim++) {
332d8f25ad8SToby Isaac     for (PetscInt deg = 1; deg <= max_deg; deg++) {
333d8f25ad8SToby Isaac       for (PetscInt form = -dim+1; form <= dim; form++) {
334*5f80ce2aSJacob Faibussowitsch         CHKERRQ(test(dim, deg, form, PetscMax(1, k), cond));
335d8f25ad8SToby Isaac       }
336d8f25ad8SToby Isaac     }
337d8f25ad8SToby Isaac   }
338d8f25ad8SToby Isaac   ierr = PetscFinalize();
339d8f25ad8SToby Isaac   return ierr;
340d8f25ad8SToby Isaac }
341d8f25ad8SToby Isaac 
342d8f25ad8SToby Isaac /*TEST
343d8f25ad8SToby Isaac 
344d8f25ad8SToby Isaac   test:
345d8f25ad8SToby Isaac     requires: !single
346d8f25ad8SToby Isaac     args:
347d8f25ad8SToby Isaac 
348d8f25ad8SToby Isaac TEST*/
349