1ccb4e88aSToby Isaac const char help[] = "Test basic creation and evaluation of PETSCSPACEPTRIMMED"; 2ccb4e88aSToby Isaac 3ccb4e88aSToby Isaac #include <petscfe.h> 4ccb4e88aSToby Isaac 5d71ae5a4SJacob Faibussowitsch static PetscErrorCode test(PetscInt dim, PetscInt formDegree, PetscInt degree, PetscInt nCopies) 6d71ae5a4SJacob Faibussowitsch { 7ccb4e88aSToby Isaac MPI_Comm comm = PETSC_COMM_SELF; 8ccb4e88aSToby Isaac PetscSpace sp; 9ccb4e88aSToby Isaac PetscInt Nf, Nb; 10ccb4e88aSToby Isaac PetscInt maxDexp, maxD, d; 11ccb4e88aSToby Isaac PetscInt Nbexp, Bsize, Dsize, Hsize; 12ccb4e88aSToby Isaac PetscReal *B, *D, *H; 13ccb4e88aSToby Isaac PetscQuadrature quad; 14ccb4e88aSToby Isaac PetscInt npoints; 15ccb4e88aSToby Isaac const PetscReal *points; 16ccb4e88aSToby Isaac 17ccb4e88aSToby Isaac PetscFunctionBegin; 189566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(comm, &sp)); 199566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)sp, "ptrimmed")); 209566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(sp, PETSCSPACEPTRIMMED)); 219566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(sp, dim)); 229566063dSJacob Faibussowitsch PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(formDegree), &Nf)); 239566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(sp, Nf * nCopies)); 249566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(sp, degree, PETSC_DETERMINE)); 259566063dSJacob Faibussowitsch PetscCall(PetscSpacePTrimmedSetFormDegree(sp, formDegree)); 269566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(sp)); 279566063dSJacob Faibussowitsch PetscCall(PetscSpaceView(sp, NULL)); 28ccb4e88aSToby Isaac 299566063dSJacob Faibussowitsch PetscCall(PetscDTPTrimmedSize(dim, formDegree == 0 ? degree : degree + 1, PetscAbsInt(formDegree), &Nbexp)); 30ccb4e88aSToby Isaac Nbexp *= nCopies; 319566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(sp, &Nb)); 3263a3b9bcSJacob Faibussowitsch PetscCheck(Nb == Nbexp, comm, PETSC_ERR_PLIB, "Space dimension mismatch, %" PetscInt_FMT " != %" PetscInt_FMT, Nbexp, Nb); 33ccb4e88aSToby Isaac 34ccb4e88aSToby Isaac maxDexp = (PetscAbsInt(formDegree) == dim || formDegree == 0) ? degree : degree + 1; 359566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(sp, &d, &maxD)); 3663a3b9bcSJacob Faibussowitsch PetscCheck(degree == d, comm, PETSC_ERR_PLIB, "Space degree mismatch, %" PetscInt_FMT " != %" PetscInt_FMT, degree, d); 3763a3b9bcSJacob Faibussowitsch PetscCheck(maxDexp == maxD, comm, PETSC_ERR_PLIB, "Space max degree mismatch, %" PetscInt_FMT " != %" PetscInt_FMT, maxDexp, maxD); 38ccb4e88aSToby Isaac 399566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, maxD + 1, -1., 1., &quad)); 409566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &npoints, &points, NULL)); 41ccb4e88aSToby Isaac 42ccb4e88aSToby Isaac Bsize = npoints * Nb * Nf * nCopies; 43ccb4e88aSToby Isaac Dsize = dim * Bsize; 44ccb4e88aSToby Isaac Hsize = dim * Dsize; 459566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(Bsize, &B, Dsize, &D, Hsize, &H)); 469566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H)); 47ad540459SPierre Jolivet for (PetscInt i = 0; i < Bsize; i++) PetscCheck(!PetscIsInfOrNanReal(B[i]), comm, PETSC_ERR_PLIB, "Bad value B[%" PetscInt_FMT "]", i); 48ad540459SPierre Jolivet for (PetscInt i = 0; i < Dsize; i++) PetscCheck(!PetscIsInfOrNanReal(D[i]), comm, PETSC_ERR_PLIB, "Bad value D[%" PetscInt_FMT "]", i); 49ad540459SPierre Jolivet for (PetscInt i = 0; i < Hsize; i++) PetscCheck(!PetscIsInfOrNanReal(H[i]), comm, PETSC_ERR_PLIB, "Bad value H[%" PetscInt_FMT "]", i); 509566063dSJacob Faibussowitsch PetscCall(PetscFree3(B, D, H)); 519566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quad)); 529566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&sp)); 53*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54ccb4e88aSToby Isaac } 55ccb4e88aSToby Isaac 56d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 57d71ae5a4SJacob Faibussowitsch { 58327415f7SBarry Smith PetscFunctionBeginUser; 599566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 60ccb4e88aSToby Isaac for (PetscInt dim = 0; dim <= 3; dim++) { 61ccb4e88aSToby Isaac for (PetscInt formDegree = -dim; formDegree <= dim; formDegree++) { 62ccb4e88aSToby Isaac for (PetscInt degree = 0; degree <= 4; degree++) { 63ccb4e88aSToby Isaac if (formDegree == 0 && degree == 0) continue; 6448a46eb9SPierre Jolivet for (PetscInt nCopies = 1; nCopies <= PetscMax(2, dim); nCopies++) PetscCall(test(dim, formDegree, degree, nCopies)); 65ccb4e88aSToby Isaac } 66ccb4e88aSToby Isaac } 67ccb4e88aSToby Isaac } 689566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 69b122ec5aSJacob Faibussowitsch return 0; 70ccb4e88aSToby Isaac } 71ccb4e88aSToby Isaac 72ccb4e88aSToby Isaac /*TEST 73ccb4e88aSToby Isaac 74ccb4e88aSToby Isaac test: 75ccb4e88aSToby Isaac 76ccb4e88aSToby Isaac TEST*/ 77