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 79371c9d4SSatish Balay 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) { 8d8f25ad8SToby Isaac PetscInt Nf; // Number of form components 9d8f25ad8SToby Isaac PetscInt Nbpt; // number of trimmed polynomials 10d8f25ad8SToby Isaac PetscInt Nk; // jet size 11d8f25ad8SToby Isaac PetscReal *p_trimmed; 12d8f25ad8SToby Isaac 13d8f25ad8SToby Isaac PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(form), &Nf)); 159566063dSJacob Faibussowitsch PetscCall(PetscDTPTrimmedSize(dim, deg, form, &Nbpt)); 169566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + jetDegree, dim, &Nk)); 179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed)); 189566063dSJacob Faibussowitsch PetscCall(PetscDTPTrimmedEvalJet(dim, npoints, points, deg, form, jetDegree, p_trimmed)); 19d8f25ad8SToby Isaac 20d8f25ad8SToby Isaac // compute the direct mass matrix 21d8f25ad8SToby Isaac PetscScalar *M_trimmed; 229566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt * Nbpt, &M_trimmed)); 23d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt; i++) { 24d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) { 25d8f25ad8SToby Isaac PetscReal v = 0.; 26d8f25ad8SToby Isaac 27d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) { 28d8f25ad8SToby Isaac const PetscReal *p_i = &p_trimmed[(i * Nf + f) * Nk * npoints]; 29d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + f) * Nk * npoints]; 30d8f25ad8SToby Isaac 319371c9d4SSatish Balay for (PetscInt pt = 0; pt < npoints; pt++) { v += p_i[pt] * p_j[pt] * weights[pt]; } 32d8f25ad8SToby Isaac } 33d8f25ad8SToby Isaac M_trimmed[i * Nbpt + j] += v; 34d8f25ad8SToby Isaac } 35d8f25ad8SToby Isaac } 36d8f25ad8SToby Isaac *_Nb = Nbpt; 37d8f25ad8SToby Isaac *_Nf = Nf; 38d8f25ad8SToby Isaac *_Nk = Nk; 39d8f25ad8SToby Isaac *B = p_trimmed; 40d8f25ad8SToby Isaac *M = M_trimmed; 41d8f25ad8SToby Isaac PetscFunctionReturn(0); 42d8f25ad8SToby Isaac } 43d8f25ad8SToby Isaac 449371c9d4SSatish Balay static PetscErrorCode test(PetscInt dim, PetscInt deg, PetscInt form, PetscInt jetDegree, PetscBool cond) { 45d8f25ad8SToby Isaac PetscQuadrature q; 46d8f25ad8SToby Isaac PetscInt npoints; 47d8f25ad8SToby Isaac const PetscReal *points; 48d8f25ad8SToby Isaac const PetscReal *weights; 49d8f25ad8SToby Isaac PetscInt Nf; // Number of form components 50d8f25ad8SToby Isaac PetscInt Nk; // jet size 51d8f25ad8SToby Isaac PetscInt Nbpt; // number of trimmed polynomials 52d8f25ad8SToby Isaac PetscReal *p_trimmed; 53d8f25ad8SToby Isaac PetscScalar *M_trimmed; 54d8f25ad8SToby Isaac PetscReal *p_scalar; 55d8f25ad8SToby Isaac PetscInt Nbp; // number of scalar polynomials 56d8f25ad8SToby Isaac PetscScalar *Mcopy; 57d8f25ad8SToby Isaac PetscScalar *M_moments; 58d8f25ad8SToby Isaac PetscReal frob_err = 0.; 59d8f25ad8SToby Isaac Mat mat_trimmed; 60d8f25ad8SToby Isaac Mat mat_moments_T; 61d8f25ad8SToby Isaac Mat AinvB; 62d8f25ad8SToby Isaac PetscInt Nbm1; 63d8f25ad8SToby Isaac Mat Mm1; 64d8f25ad8SToby Isaac PetscReal *p_trimmed_copy; 65d8f25ad8SToby Isaac PetscReal *M_moment_real; 66d8f25ad8SToby Isaac 67d8f25ad8SToby Isaac PetscFunctionBegin; 68d8f25ad8SToby Isaac // Construct an appropriate quadrature 699566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, deg + 2, -1., 1., &q)); 709566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights)); 71d8f25ad8SToby Isaac 729566063dSJacob Faibussowitsch PetscCall(constructTabulationAndMass(dim, deg, form, jetDegree, npoints, points, weights, &Nbpt, &Nf, &Nk, &p_trimmed, &M_trimmed)); 73d8f25ad8SToby Isaac 749566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + deg, dim, &Nbp)); 759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbp * Nk * npoints, &p_scalar)); 769566063dSJacob Faibussowitsch PetscCall(PetscDTPKDEvalJet(dim, npoints, points, deg, jetDegree, p_scalar)); 77d8f25ad8SToby Isaac 789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt * Nbpt, &Mcopy)); 79d8f25ad8SToby Isaac // Print the condition numbers (useful for testing out different bases internally in PetscDTPTrimmedEvalJet()) 80d8f25ad8SToby Isaac #if !defined(PETSC_USE_COMPLEX) 81d8f25ad8SToby Isaac if (cond) { 82d8f25ad8SToby Isaac PetscReal *S; 83d8f25ad8SToby Isaac PetscScalar *work; 84d8f25ad8SToby Isaac PetscBLASInt n = Nbpt; 85d8f25ad8SToby Isaac PetscBLASInt lwork = 5 * Nbpt; 86d8f25ad8SToby Isaac PetscBLASInt lierr; 87d8f25ad8SToby Isaac 889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt, &S)); 899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(5 * Nbpt, &work)); 909566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt)); 91d8f25ad8SToby Isaac 92792fecdfSBarry Smith PetscCallBLAS("LAPACKgesvd", LAPACKgesvd_("N", "N", &n, &n, Mcopy, &n, S, NULL, &n, NULL, &n, work, &lwork, &lierr)); 93d8f25ad8SToby Isaac PetscReal cond = S[0] / S[Nbpt - 1]; 9463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "dimension %" PetscInt_FMT ", degree %" PetscInt_FMT ", form %" PetscInt_FMT ": condition number %g\n", dim, deg, form, (double)cond)); 959566063dSJacob Faibussowitsch PetscCall(PetscFree(work)); 969566063dSJacob Faibussowitsch PetscCall(PetscFree(S)); 97d8f25ad8SToby Isaac } 98d8f25ad8SToby Isaac #endif 99d8f25ad8SToby Isaac 100d8f25ad8SToby Isaac // compute the moments with the orthonormal polynomials 1019566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt * Nbp * Nf, &M_moments)); 102d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbp; i++) { 103d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) { 104d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) { 105d8f25ad8SToby Isaac PetscReal v = 0.; 106d8f25ad8SToby Isaac const PetscReal *p_i = &p_scalar[i * Nk * npoints]; 107d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + f) * Nk * npoints]; 108d8f25ad8SToby Isaac 1099371c9d4SSatish Balay for (PetscInt pt = 0; pt < npoints; pt++) { v += p_i[pt] * p_j[pt] * weights[pt]; } 110d8f25ad8SToby Isaac M_moments[(i * Nf + f) * Nbpt + j] += v; 111d8f25ad8SToby Isaac } 112d8f25ad8SToby Isaac } 113d8f25ad8SToby Isaac } 114d8f25ad8SToby Isaac 115d8f25ad8SToby Isaac // subtract M_moments^T * M_moments from M_trimmed: because the trimmed polynomials should be contained in 116d8f25ad8SToby Isaac // the full polynomials, the result should be zero 1179566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt)); 118d8f25ad8SToby Isaac { 119d8f25ad8SToby Isaac PetscBLASInt m = Nbpt; 120d8f25ad8SToby Isaac PetscBLASInt n = Nbpt; 121d8f25ad8SToby Isaac PetscBLASInt k = Nbp * Nf; 122d8f25ad8SToby Isaac PetscScalar mone = -1.; 123d8f25ad8SToby Isaac PetscScalar one = 1.; 124d8f25ad8SToby Isaac 125792fecdfSBarry Smith PetscCallBLAS("BLASgemm", BLASgemm_("N", "T", &m, &n, &k, &mone, M_moments, &m, M_moments, &m, &one, Mcopy, &m)); 126d8f25ad8SToby Isaac } 127d8f25ad8SToby Isaac 128d8f25ad8SToby Isaac frob_err = 0.; 129d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt * Nbpt; i++) frob_err += PetscRealPart(Mcopy[i]) * PetscRealPart(Mcopy[i]); 130d8f25ad8SToby Isaac frob_err = PetscSqrtReal(frob_err); 131d8f25ad8SToby Isaac 1327a46b595SBarry 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); 133d8f25ad8SToby Isaac 134d8f25ad8SToby 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 1359566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt, M_trimmed, &mat_trimmed)); 1369566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim + deg - 1, dim, &Nbm1)); 1379566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbm1 * Nf, M_moments, &mat_moments_T)); 1389566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat_moments_T, MAT_DO_NOT_COPY_VALUES, &AinvB)); 1399566063dSJacob Faibussowitsch PetscCall(MatLUFactor(mat_trimmed, NULL, NULL, NULL)); 1409566063dSJacob Faibussowitsch PetscCall(MatMatSolve(mat_trimmed, mat_moments_T, AinvB)); 1419566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(mat_moments_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Mm1)); 1429566063dSJacob Faibussowitsch PetscCall(MatShift(Mm1, -1.)); 1439566063dSJacob Faibussowitsch PetscCall(MatNorm(Mm1, NORM_FROBENIUS, &frob_err)); 1447a46b595SBarry 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); 1459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Mm1)); 1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AinvB)); 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_moments_T)); 148d8f25ad8SToby Isaac 149d8f25ad8SToby Isaac // The Koszul differential applied to P trimmed (Lambda k+1) should be contained in P trimmed (Lambda k) 150d8f25ad8SToby Isaac if (PetscAbsInt(form) < dim) { 151d8f25ad8SToby Isaac PetscInt Nf1, Nbpt1, Nk1; 152d8f25ad8SToby Isaac PetscReal *p_trimmed1; 153d8f25ad8SToby Isaac PetscScalar *M_trimmed1; 154d8f25ad8SToby Isaac PetscInt(*pattern)[3]; 155d8f25ad8SToby Isaac PetscReal *p_koszul; 156d8f25ad8SToby Isaac PetscScalar *M_koszul; 157d8f25ad8SToby Isaac PetscScalar *M_k_moment; 158d8f25ad8SToby Isaac Mat mat_koszul; 159d8f25ad8SToby Isaac Mat mat_k_moment_T; 160d8f25ad8SToby Isaac Mat AinvB; 161d8f25ad8SToby Isaac Mat prod; 162d8f25ad8SToby Isaac 1639371c9d4SSatish Balay PetscCall(constructTabulationAndMass(dim, deg, form < 0 ? form - 1 : form + 1, 0, npoints, points, weights, &Nbpt1, &Nf1, &Nk1, &p_trimmed1, &M_trimmed1)); 164d8f25ad8SToby Isaac 1659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nf1 * (PetscAbsInt(form) + 1), &pattern)); 1669566063dSJacob Faibussowitsch PetscCall(PetscDTAltVInteriorPattern(dim, PetscAbsInt(form) + 1, pattern)); 167d8f25ad8SToby Isaac 168d8f25ad8SToby Isaac // apply the Koszul operator 1699566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt1 * Nf * npoints, &p_koszul)); 170d8f25ad8SToby Isaac for (PetscInt b = 0; b < Nbpt1; b++) { 171d8f25ad8SToby Isaac for (PetscInt a = 0; a < Nf1 * (PetscAbsInt(form) + 1); a++) { 172d8f25ad8SToby Isaac PetscInt i, j, k; 173d8f25ad8SToby Isaac PetscReal sign; 174d8f25ad8SToby Isaac PetscReal *p_i; 175d8f25ad8SToby Isaac const PetscReal *p_j; 176d8f25ad8SToby Isaac 177d8f25ad8SToby Isaac i = pattern[a][0]; 1789371c9d4SSatish Balay if (form < 0) { i = Nf - 1 - i; } 179d8f25ad8SToby Isaac j = pattern[a][1]; 1809371c9d4SSatish Balay if (form < 0) { j = Nf1 - 1 - j; } 181d8f25ad8SToby Isaac k = pattern[a][2] < 0 ? -(pattern[a][2] + 1) : pattern[a][2]; 182d8f25ad8SToby Isaac sign = pattern[a][2] < 0 ? -1 : 1; 1839371c9d4SSatish Balay if (form < 0 && (i & 1) ^ (j & 1)) { sign = -sign; } 184d8f25ad8SToby Isaac 185d8f25ad8SToby Isaac p_i = &p_koszul[(b * Nf + i) * npoints]; 186d8f25ad8SToby Isaac p_j = &p_trimmed1[(b * Nf1 + j) * npoints]; 1879371c9d4SSatish Balay for (PetscInt pt = 0; pt < npoints; pt++) { p_i[pt] += p_j[pt] * points[pt * dim + k] * sign; } 188d8f25ad8SToby Isaac } 189d8f25ad8SToby Isaac } 190d8f25ad8SToby Isaac 191d8f25ad8SToby Isaac // mass matrix of the result 1929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt1 * Nbpt1, &M_koszul)); 193d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt1; i++) { 194d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt1; j++) { 195d8f25ad8SToby Isaac PetscReal val = 0.; 196d8f25ad8SToby Isaac 197d8f25ad8SToby Isaac for (PetscInt v = 0; v < Nf; v++) { 198d8f25ad8SToby Isaac const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints]; 199d8f25ad8SToby Isaac const PetscReal *p_j = &p_koszul[(j * Nf + v) * npoints]; 200d8f25ad8SToby Isaac 2019371c9d4SSatish Balay for (PetscInt pt = 0; pt < npoints; pt++) { val += p_i[pt] * p_j[pt] * weights[pt]; } 202d8f25ad8SToby Isaac } 203d8f25ad8SToby Isaac M_koszul[i * Nbpt1 + j] = val; 204d8f25ad8SToby Isaac } 205d8f25ad8SToby Isaac } 206d8f25ad8SToby Isaac 207d8f25ad8SToby Isaac // moment matrix between the result and P trimmed 2089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbpt * Nbpt1, &M_k_moment)); 209d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt1; i++) { 210d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) { 211d8f25ad8SToby Isaac PetscReal val = 0.; 212d8f25ad8SToby Isaac 213d8f25ad8SToby Isaac for (PetscInt v = 0; v < Nf; v++) { 214d8f25ad8SToby Isaac const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints]; 215d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + v) * Nk * npoints]; 216d8f25ad8SToby Isaac 2179371c9d4SSatish Balay for (PetscInt pt = 0; pt < npoints; pt++) { val += p_i[pt] * p_j[pt] * weights[pt]; } 218d8f25ad8SToby Isaac } 219d8f25ad8SToby Isaac M_k_moment[i * Nbpt + j] = val; 220d8f25ad8SToby Isaac } 221d8f25ad8SToby Isaac } 222d8f25ad8SToby Isaac 223d8f25ad8SToby Isaac // M_k_moment M_trimmed^{-1} M_k_moment^T == M_koszul 2249566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt1, Nbpt1, M_koszul, &mat_koszul)); 2259566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt1, M_k_moment, &mat_k_moment_T)); 2269566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mat_k_moment_T, MAT_DO_NOT_COPY_VALUES, &AinvB)); 2279566063dSJacob Faibussowitsch PetscCall(MatMatSolve(mat_trimmed, mat_k_moment_T, AinvB)); 2289566063dSJacob Faibussowitsch PetscCall(MatTransposeMatMult(mat_k_moment_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &prod)); 2299566063dSJacob Faibussowitsch PetscCall(MatAXPY(prod, -1., mat_koszul, SAME_NONZERO_PATTERN)); 2309566063dSJacob Faibussowitsch PetscCall(MatNorm(prod, NORM_FROBENIUS, &frob_err)); 231d8f25ad8SToby Isaac if (frob_err > PETSC_SMALL) { 23263a3b9bcSJacob Faibussowitsch SETERRQ(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); 233d8f25ad8SToby Isaac } 234d8f25ad8SToby Isaac 2359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&prod)); 2369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AinvB)); 2379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_k_moment_T)); 2389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_koszul)); 2399566063dSJacob Faibussowitsch PetscCall(PetscFree(M_k_moment)); 2409566063dSJacob Faibussowitsch PetscCall(PetscFree(M_koszul)); 2419566063dSJacob Faibussowitsch PetscCall(PetscFree(p_koszul)); 2429566063dSJacob Faibussowitsch PetscCall(PetscFree(pattern)); 2439566063dSJacob Faibussowitsch PetscCall(PetscFree(p_trimmed1)); 2449566063dSJacob Faibussowitsch PetscCall(PetscFree(M_trimmed1)); 245d8f25ad8SToby Isaac } 246d8f25ad8SToby Isaac 247d8f25ad8SToby Isaac // M_moments has shape [Nbp][Nf][Nbpt] 248d8f25ad8SToby Isaac // p_scalar has shape [Nbp][Nk][npoints] 249d8f25ad8SToby Isaac // contracting on [Nbp] should be the same shape as 250d8f25ad8SToby Isaac // p_trimmed, which is [Nbpt][Nf][Nk][npoints] 2519566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed_copy)); 2529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nbp * Nf * Nbpt, &M_moment_real)); 2539371c9d4SSatish Balay for (PetscInt i = 0; i < Nbp * Nf * Nbpt; i++) { M_moment_real[i] = PetscRealPart(M_moments[i]); } 254d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) { 255d8f25ad8SToby Isaac PetscBLASInt m = Nk * npoints; 256d8f25ad8SToby Isaac PetscBLASInt n = Nbpt; 257d8f25ad8SToby Isaac PetscBLASInt k = Nbp; 258d8f25ad8SToby Isaac PetscBLASInt lda = Nk * npoints; 259d8f25ad8SToby Isaac PetscBLASInt ldb = Nf * Nbpt; 260d8f25ad8SToby Isaac PetscBLASInt ldc = Nf * Nk * npoints; 261d8f25ad8SToby Isaac PetscReal alpha = 1.0; 262d8f25ad8SToby Isaac PetscReal beta = 1.0; 263d8f25ad8SToby Isaac 264792fecdfSBarry 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)); 265d8f25ad8SToby Isaac } 266d8f25ad8SToby Isaac frob_err = 0.; 2679371c9d4SSatish Balay 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]); } 268d8f25ad8SToby Isaac frob_err = PetscSqrtReal(frob_err); 269d8f25ad8SToby Isaac 2701baa6e33SBarry 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); 271d8f25ad8SToby Isaac 2729566063dSJacob Faibussowitsch PetscCall(PetscFree(M_moment_real)); 2739566063dSJacob Faibussowitsch PetscCall(PetscFree(p_trimmed_copy)); 2749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat_trimmed)); 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(Mcopy)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree(M_moments)); 2779566063dSJacob Faibussowitsch PetscCall(PetscFree(M_trimmed)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFree(p_trimmed)); 2799566063dSJacob Faibussowitsch PetscCall(PetscFree(p_scalar)); 2809566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 281d8f25ad8SToby Isaac PetscFunctionReturn(0); 282d8f25ad8SToby Isaac } 283d8f25ad8SToby Isaac 2849371c9d4SSatish Balay int main(int argc, char **argv) { 285d8f25ad8SToby Isaac PetscInt max_dim = 3; 286d8f25ad8SToby Isaac PetscInt max_deg = 4; 287d8f25ad8SToby Isaac PetscInt k = 3; 288d8f25ad8SToby Isaac PetscBool cond = PETSC_FALSE; 289d8f25ad8SToby Isaac 290327415f7SBarry Smith PetscFunctionBeginUser; 2919566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 292d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for PetscDTPTrimmedEvalJet() tests", "none"); 2939566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_dim", "Maximum dimension of the simplex", __FILE__, max_dim, &max_dim, NULL)); 2949566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_degree", "Maximum degree of the trimmed polynomial space", __FILE__, max_deg, &max_deg, NULL)); 2959566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_jet", "The number of derivatives to test", __FILE__, k, &k, NULL)); 2969566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-cond", "Compute the condition numbers of the mass matrices of the bases", __FILE__, cond, &cond, NULL)); 297d0609cedSBarry Smith PetscOptionsEnd(); 298d8f25ad8SToby Isaac for (PetscInt dim = 2; dim <= max_dim; dim++) { 299d8f25ad8SToby Isaac for (PetscInt deg = 1; deg <= max_deg; deg++) { 300*48a46eb9SPierre Jolivet for (PetscInt form = -dim + 1; form <= dim; form++) PetscCall(test(dim, deg, form, PetscMax(1, k), cond)); 301d8f25ad8SToby Isaac } 302d8f25ad8SToby Isaac } 3039566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 304b122ec5aSJacob Faibussowitsch return 0; 305d8f25ad8SToby Isaac } 306d8f25ad8SToby Isaac 307d8f25ad8SToby Isaac /*TEST 308d8f25ad8SToby Isaac 309d8f25ad8SToby Isaac test: 310d8f25ad8SToby Isaac requires: !single 311d8f25ad8SToby Isaac args: 312d8f25ad8SToby Isaac 313d8f25ad8SToby Isaac TEST*/ 314