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