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; 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(dim, PetscAbsInt(form), &Nf)); 195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTPTrimmedSize(dim, deg, form, &Nbpt)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(dim + jetDegree, dim, &Nk)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed)); 225f80ce2aSJacob 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; 265f80ce2aSJacob 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 74d8f25ad8SToby Isaac PetscFunctionBegin; 75d8f25ad8SToby Isaac // Construct an appropriate quadrature 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTStroudConicalQuadrature(dim, 1, deg + 2, -1., 1., &q)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureGetData(q, NULL, NULL, &npoints, &points, &weights)); 78d8f25ad8SToby Isaac 795f80ce2aSJacob Faibussowitsch CHKERRQ(constructTabulationAndMass(dim, deg, form, jetDegree, npoints, points, weights, &Nbpt, &Nf, &Nk, &p_trimmed, &M_trimmed)); 80d8f25ad8SToby Isaac 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(dim + deg, dim, &Nbp)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nbp * Nk * npoints, &p_scalar)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTPKDEvalJet(dim, npoints, points, deg, jetDegree, p_scalar)); 84d8f25ad8SToby Isaac 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nbpt * Nbpt, &Mcopy)); 86d8f25ad8SToby Isaac // Print the condition numbers (useful for testing out different bases internally in PetscDTPTrimmedEvalJet()) 87d8f25ad8SToby Isaac #if !defined(PETSC_USE_COMPLEX) 88d8f25ad8SToby Isaac if (cond) { 89d8f25ad8SToby Isaac PetscReal *S; 90d8f25ad8SToby Isaac PetscScalar *work; 91d8f25ad8SToby Isaac PetscBLASInt n = Nbpt; 92d8f25ad8SToby Isaac PetscBLASInt lwork = 5 * Nbpt; 93d8f25ad8SToby Isaac PetscBLASInt lierr; 94d8f25ad8SToby Isaac 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nbpt, &S)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(5*Nbpt, &work)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt)); 98d8f25ad8SToby Isaac 99d8f25ad8SToby Isaac PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("N","N",&n,&n,Mcopy,&n,S,NULL,&n,NULL,&n,work,&lwork,&lierr)); 100d8f25ad8SToby Isaac PetscReal cond = S[0] / S[Nbpt - 1]; 101*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "dimension %D, degree %D, form %D: condition number %g\n", dim, deg, form, (double) cond)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(work)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(S)); 104d8f25ad8SToby Isaac } 105d8f25ad8SToby Isaac #endif 106d8f25ad8SToby Isaac 107d8f25ad8SToby Isaac // compute the moments with the orthonormal polynomials 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(Nbpt * Nbp * Nf, &M_moments)); 109d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbp; i++) { 110d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) { 111d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) { 112d8f25ad8SToby Isaac PetscReal v = 0.; 113d8f25ad8SToby Isaac const PetscReal *p_i = &p_scalar[i * Nk * npoints]; 114d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + f) * Nk * npoints]; 115d8f25ad8SToby Isaac 116d8f25ad8SToby Isaac for (PetscInt pt = 0; pt < npoints; pt++) { 117d8f25ad8SToby Isaac v += p_i[pt] * p_j[pt] * weights[pt]; 118d8f25ad8SToby Isaac } 119d8f25ad8SToby Isaac M_moments[(i * Nf + f) * Nbpt + j] += v; 120d8f25ad8SToby Isaac } 121d8f25ad8SToby Isaac } 122d8f25ad8SToby Isaac } 123d8f25ad8SToby Isaac 124d8f25ad8SToby Isaac // subtract M_moments^T * M_moments from M_trimmed: because the trimmed polynomials should be contained in 125d8f25ad8SToby Isaac // the full polynomials, the result should be zero 1265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(Mcopy, M_trimmed, Nbpt * Nbpt)); 127d8f25ad8SToby Isaac { 128d8f25ad8SToby Isaac PetscBLASInt m = Nbpt; 129d8f25ad8SToby Isaac PetscBLASInt n = Nbpt; 130d8f25ad8SToby Isaac PetscBLASInt k = Nbp * Nf; 131d8f25ad8SToby Isaac PetscScalar mone = -1.; 132d8f25ad8SToby Isaac PetscScalar one = 1.; 133d8f25ad8SToby Isaac 134d8f25ad8SToby Isaac PetscStackCallBLAS("BLASgemm",BLASgemm_("N","T",&m,&n,&k,&mone,M_moments,&m,M_moments,&m,&one,Mcopy,&m)); 135d8f25ad8SToby Isaac } 136d8f25ad8SToby Isaac 137d8f25ad8SToby Isaac frob_err = 0.; 138d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt * Nbpt; i++) frob_err += PetscRealPart(Mcopy[i]) * PetscRealPart(Mcopy[i]); 139d8f25ad8SToby Isaac frob_err = PetscSqrtReal(frob_err); 140d8f25ad8SToby Isaac 141d8f25ad8SToby Isaac if (frob_err > PETSC_SMALL) { 14298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %D, degree %D, form %D: trimmed projection error %g", dim, deg, form, (double) frob_err); 143d8f25ad8SToby Isaac } 144d8f25ad8SToby Isaac 145d8f25ad8SToby 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 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt, M_trimmed, &mat_trimmed)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTBinomialInt(dim + deg - 1, dim, &Nbm1)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbm1 * Nf, M_moments, &mat_moments_T)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(mat_moments_T, MAT_DO_NOT_COPY_VALUES, &AinvB)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactor(mat_trimmed, NULL, NULL, NULL)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(mat_trimmed, mat_moments_T, AinvB)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(mat_moments_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &Mm1)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(Mm1, -1.)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(Mm1, NORM_FROBENIUS, &frob_err)); 155d8f25ad8SToby Isaac if (frob_err > PETSC_SMALL) { 15698921bdaSJacob 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); 157d8f25ad8SToby Isaac } 1585f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Mm1)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&AinvB)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat_moments_T)); 161d8f25ad8SToby Isaac 162d8f25ad8SToby Isaac // The Koszul differential applied to P trimmed (Lambda k+1) should be contained in P trimmed (Lambda k) 163d8f25ad8SToby Isaac if (PetscAbsInt(form) < dim) { 164d8f25ad8SToby Isaac PetscInt Nf1, Nbpt1, Nk1; 165d8f25ad8SToby Isaac PetscReal *p_trimmed1; 166d8f25ad8SToby Isaac PetscScalar *M_trimmed1; 167d8f25ad8SToby Isaac PetscInt (*pattern)[3]; 168d8f25ad8SToby Isaac PetscReal *p_koszul; 169d8f25ad8SToby Isaac PetscScalar *M_koszul; 170d8f25ad8SToby Isaac PetscScalar *M_k_moment; 171d8f25ad8SToby Isaac Mat mat_koszul; 172d8f25ad8SToby Isaac Mat mat_k_moment_T; 173d8f25ad8SToby Isaac Mat AinvB; 174d8f25ad8SToby Isaac Mat prod; 175d8f25ad8SToby Isaac 176*b122ec5aSJacob Faibussowitsch CHKERRQ(constructTabulationAndMass(dim, deg, form < 0 ? form - 1 : form + 1, 0, npoints, points, weights, &Nbpt1, &Nf1, &Nk1, 177*b122ec5aSJacob Faibussowitsch &p_trimmed1, &M_trimmed1)); 178d8f25ad8SToby Isaac 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nf1 * (PetscAbsInt(form) + 1), &pattern)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscDTAltVInteriorPattern(dim, PetscAbsInt(form) + 1, pattern)); 181d8f25ad8SToby Isaac 182d8f25ad8SToby Isaac // apply the Koszul operator 1835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(Nbpt1 * Nf * npoints, &p_koszul)); 184d8f25ad8SToby Isaac for (PetscInt b = 0; b < Nbpt1; b++) { 185d8f25ad8SToby Isaac for (PetscInt a = 0; a < Nf1 * (PetscAbsInt(form) + 1); a++) { 186d8f25ad8SToby Isaac PetscInt i,j,k; 187d8f25ad8SToby Isaac PetscReal sign; 188d8f25ad8SToby Isaac PetscReal *p_i; 189d8f25ad8SToby Isaac const PetscReal *p_j; 190d8f25ad8SToby Isaac 191d8f25ad8SToby Isaac i = pattern[a][0]; 192d8f25ad8SToby Isaac if (form < 0) { 193d8f25ad8SToby Isaac i = Nf-1-i; 194d8f25ad8SToby Isaac } 195d8f25ad8SToby Isaac j = pattern[a][1]; 196d8f25ad8SToby Isaac if (form < 0) { 197d8f25ad8SToby Isaac j = Nf1-1-j; 198d8f25ad8SToby Isaac } 199d8f25ad8SToby Isaac k = pattern[a][2] < 0 ? -(pattern[a][2] + 1) : pattern[a][2]; 200d8f25ad8SToby Isaac sign = pattern[a][2] < 0 ? -1 : 1; 201d8f25ad8SToby Isaac if (form < 0 && (i & 1) ^ (j & 1)) { 202d8f25ad8SToby Isaac sign = -sign; 203d8f25ad8SToby Isaac } 204d8f25ad8SToby Isaac 205d8f25ad8SToby Isaac p_i = &p_koszul[(b * Nf + i) * npoints]; 206d8f25ad8SToby Isaac p_j = &p_trimmed1[(b * Nf1 + j) * npoints]; 207d8f25ad8SToby Isaac for (PetscInt pt = 0; pt < npoints; pt++) { 208d8f25ad8SToby Isaac p_i[pt] += p_j[pt] * points[pt * dim + k] * sign; 209d8f25ad8SToby Isaac } 210d8f25ad8SToby Isaac } 211d8f25ad8SToby Isaac } 212d8f25ad8SToby Isaac 213d8f25ad8SToby Isaac // mass matrix of the result 2145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nbpt1 * Nbpt1, &M_koszul)); 215d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt1; i++) { 216d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt1; j++) { 217d8f25ad8SToby Isaac PetscReal val = 0.; 218d8f25ad8SToby Isaac 219d8f25ad8SToby Isaac for (PetscInt v = 0; v < Nf; v++) { 220d8f25ad8SToby Isaac const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints]; 221d8f25ad8SToby Isaac const PetscReal *p_j = &p_koszul[(j * Nf + v) * npoints]; 222d8f25ad8SToby Isaac 223d8f25ad8SToby Isaac for (PetscInt pt = 0; pt < npoints; pt++) { 224d8f25ad8SToby Isaac val += p_i[pt] * p_j[pt] * weights[pt]; 225d8f25ad8SToby Isaac } 226d8f25ad8SToby Isaac } 227d8f25ad8SToby Isaac M_koszul[i * Nbpt1 + j] = val; 228d8f25ad8SToby Isaac } 229d8f25ad8SToby Isaac } 230d8f25ad8SToby Isaac 231d8f25ad8SToby Isaac // moment matrix between the result and P trimmed 2325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nbpt * Nbpt1, &M_k_moment)); 233d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt1; i++) { 234d8f25ad8SToby Isaac for (PetscInt j = 0; j < Nbpt; j++) { 235d8f25ad8SToby Isaac PetscReal val = 0.; 236d8f25ad8SToby Isaac 237d8f25ad8SToby Isaac for (PetscInt v = 0; v < Nf; v++) { 238d8f25ad8SToby Isaac const PetscReal *p_i = &p_koszul[(i * Nf + v) * npoints]; 239d8f25ad8SToby Isaac const PetscReal *p_j = &p_trimmed[(j * Nf + v) * Nk * npoints]; 240d8f25ad8SToby Isaac 241d8f25ad8SToby Isaac for (PetscInt pt = 0; pt < npoints; pt++) { 242d8f25ad8SToby Isaac val += p_i[pt] * p_j[pt] * weights[pt]; 243d8f25ad8SToby Isaac } 244d8f25ad8SToby Isaac } 245d8f25ad8SToby Isaac M_k_moment[i * Nbpt + j] = val; 246d8f25ad8SToby Isaac } 247d8f25ad8SToby Isaac } 248d8f25ad8SToby Isaac 249d8f25ad8SToby Isaac // M_k_moment M_trimmed^{-1} M_k_moment^T == M_koszul 2505f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt1, Nbpt1, M_koszul, &mat_koszul)); 2515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSeqDense(PETSC_COMM_SELF, Nbpt, Nbpt1, M_k_moment, &mat_k_moment_T)); 2525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(mat_k_moment_T, MAT_DO_NOT_COPY_VALUES, &AinvB)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(mat_trimmed, mat_k_moment_T, AinvB)); 2545f80ce2aSJacob Faibussowitsch CHKERRQ(MatTransposeMatMult(mat_k_moment_T, AinvB, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &prod)); 2555f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(prod, -1., mat_koszul, SAME_NONZERO_PATTERN)); 2565f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(prod, NORM_FROBENIUS, &frob_err)); 257d8f25ad8SToby Isaac if (frob_err > PETSC_SMALL) { 25898921bdaSJacob 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); 259d8f25ad8SToby Isaac } 260d8f25ad8SToby Isaac 2615f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&prod)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&AinvB)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat_k_moment_T)); 2645f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat_koszul)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(M_k_moment)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(M_koszul)); 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(p_koszul)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pattern)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(p_trimmed1)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(M_trimmed1)); 271d8f25ad8SToby Isaac } 272d8f25ad8SToby Isaac 273d8f25ad8SToby Isaac // M_moments has shape [Nbp][Nf][Nbpt] 274d8f25ad8SToby Isaac // p_scalar has shape [Nbp][Nk][npoints] 275d8f25ad8SToby Isaac // contracting on [Nbp] should be the same shape as 276d8f25ad8SToby Isaac // p_trimmed, which is [Nbpt][Nf][Nk][npoints] 2775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(Nbpt * Nf * Nk * npoints, &p_trimmed_copy)); 2785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(Nbp * Nf * Nbpt, &M_moment_real)); 279d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbp * Nf * Nbpt; i++) { 280d8f25ad8SToby Isaac M_moment_real[i] = PetscRealPart(M_moments[i]); 281d8f25ad8SToby Isaac } 282d8f25ad8SToby Isaac for (PetscInt f = 0; f < Nf; f++) { 283d8f25ad8SToby Isaac PetscBLASInt m = Nk * npoints; 284d8f25ad8SToby Isaac PetscBLASInt n = Nbpt; 285d8f25ad8SToby Isaac PetscBLASInt k = Nbp; 286d8f25ad8SToby Isaac PetscBLASInt lda = Nk * npoints; 287d8f25ad8SToby Isaac PetscBLASInt ldb = Nf * Nbpt; 288d8f25ad8SToby Isaac PetscBLASInt ldc = Nf * Nk * npoints; 289d8f25ad8SToby Isaac PetscReal alpha = 1.0; 290d8f25ad8SToby Isaac PetscReal beta = 1.0; 291d8f25ad8SToby Isaac 292d8f25ad8SToby 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)); 293d8f25ad8SToby Isaac } 294d8f25ad8SToby Isaac frob_err = 0.; 295d8f25ad8SToby Isaac for (PetscInt i = 0; i < Nbpt * Nf * Nk * npoints; i++) { 296d8f25ad8SToby Isaac frob_err += (p_trimmed_copy[i] - p_trimmed[i]) * (p_trimmed_copy[i] - p_trimmed[i]); 297d8f25ad8SToby Isaac } 298d8f25ad8SToby Isaac frob_err = PetscSqrtReal(frob_err); 299d8f25ad8SToby Isaac 300d8f25ad8SToby Isaac if (frob_err > PETSC_SMALL) { 30198921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_PLIB, "dimension %D, degree %D, form %D: jet error %g", dim, deg, form, (double) frob_err); 302d8f25ad8SToby Isaac } 303d8f25ad8SToby Isaac 3045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(M_moment_real)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(p_trimmed_copy)); 3065f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat_trimmed)); 3075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(Mcopy)); 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(M_moments)); 3095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(M_trimmed)); 3105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(p_trimmed)); 3115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(p_scalar)); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscQuadratureDestroy(&q)); 313d8f25ad8SToby Isaac PetscFunctionReturn(0); 314d8f25ad8SToby Isaac } 315d8f25ad8SToby Isaac 316d8f25ad8SToby Isaac int main(int argc, char **argv) 317d8f25ad8SToby Isaac { 318d8f25ad8SToby Isaac PetscInt max_dim = 3; 319d8f25ad8SToby Isaac PetscInt max_deg = 4; 320d8f25ad8SToby Isaac PetscInt k = 3; 321d8f25ad8SToby Isaac PetscBool cond = PETSC_FALSE; 322*b122ec5aSJacob Faibussowitsch PetscErrorCode ierr; 323d8f25ad8SToby Isaac 324*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc, &argv, NULL, help)); 325d8f25ad8SToby Isaac ierr = PetscOptionsBegin(PETSC_COMM_WORLD,"","Options for PetscDTPTrimmedEvalJet() tests","none");CHKERRQ(ierr); 3265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-max_dim", "Maximum dimension of the simplex",__FILE__,max_dim,&max_dim,NULL)); 3275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-max_degree", "Maximum degree of the trimmed polynomial space",__FILE__,max_deg,&max_deg,NULL)); 3285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-max_jet", "The number of derivatives to test",__FILE__,k,&k,NULL)); 3295f80ce2aSJacob 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++) { 3345f80ce2aSJacob Faibussowitsch CHKERRQ(test(dim, deg, form, PetscMax(1, k), cond)); 335d8f25ad8SToby Isaac } 336d8f25ad8SToby Isaac } 337d8f25ad8SToby Isaac } 338*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 339*b122ec5aSJacob Faibussowitsch return 0; 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