11bd0d7efSStefano Zampini static const char help[] = "Tests DMCreateMassMatrix and DMCreateMassMatrixLumped"; 21bd0d7efSStefano Zampini 31bd0d7efSStefano Zampini #include <petscdmplex.h> 41bd0d7efSStefano Zampini #include <petscfe.h> 51bd0d7efSStefano Zampini #include <petscds.h> 61bd0d7efSStefano Zampini 71bd0d7efSStefano Zampini static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[]) 81bd0d7efSStefano Zampini { 91bd0d7efSStefano Zampini obj[0] = 1.0; 101bd0d7efSStefano Zampini } 111bd0d7efSStefano Zampini 121bd0d7efSStefano Zampini static PetscErrorCode SetupDiscretization(DM dm) 131bd0d7efSStefano Zampini { 141bd0d7efSStefano Zampini PetscFE fe0, fe1, fe2; 151bd0d7efSStefano Zampini DM plex; 161bd0d7efSStefano Zampini PetscInt dim; 171bd0d7efSStefano Zampini PetscBool simplex; 181bd0d7efSStefano Zampini PetscDS ds; 191bd0d7efSStefano Zampini 201bd0d7efSStefano Zampini PetscFunctionBeginUser; 211bd0d7efSStefano Zampini PetscCall(DMGetDimension(dm, &dim)); 221bd0d7efSStefano Zampini PetscCall(DMConvert(dm, DMPLEX, &plex)); 231bd0d7efSStefano Zampini PetscCall(DMPlexIsSimplex(plex, &simplex)); 24*5440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &simplex, 1, MPI_C_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm))); 251bd0d7efSStefano Zampini PetscCall(DMDestroy(&plex)); 261bd0d7efSStefano Zampini 271bd0d7efSStefano Zampini /* Create finite element */ 281bd0d7efSStefano Zampini PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, "f0_", -1, &fe0)); 291bd0d7efSStefano Zampini PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 2, simplex, "f1_", -1, &fe1)); 301bd0d7efSStefano Zampini PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, simplex, "f2_", -1, &fe2)); 311bd0d7efSStefano Zampini PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe0)); 321bd0d7efSStefano Zampini PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe1)); 331bd0d7efSStefano Zampini PetscCall(DMSetField(dm, 2, NULL, (PetscObject)fe2)); 341bd0d7efSStefano Zampini PetscCall(DMCreateDS(dm)); 351bd0d7efSStefano Zampini PetscCall(DMGetDS(dm, &ds)); 361bd0d7efSStefano Zampini PetscCall(PetscDSSetObjective(ds, 0, volume)); 371bd0d7efSStefano Zampini PetscCall(PetscDSSetObjective(ds, 1, volume)); 381bd0d7efSStefano Zampini PetscCall(PetscDSSetObjective(ds, 2, volume)); 391bd0d7efSStefano Zampini PetscCall(PetscFEDestroy(&fe0)); 401bd0d7efSStefano Zampini PetscCall(PetscFEDestroy(&fe1)); 411bd0d7efSStefano Zampini PetscCall(PetscFEDestroy(&fe2)); 421bd0d7efSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 431bd0d7efSStefano Zampini } 441bd0d7efSStefano Zampini 451bd0d7efSStefano Zampini #define CheckVals(a, b, rtol, atol, msg) \ 461bd0d7efSStefano Zampini do { \ 471bd0d7efSStefano Zampini if (!PetscIsCloseAtTolScalar(a, b, rtol, atol)) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s: %g (%g - %g)\n", msg, (double)PetscAbsScalar(a - b), (double)PetscAbsScalar(a), (double)PetscAbsScalar(b))); \ 481bd0d7efSStefano Zampini } while (0) 491bd0d7efSStefano Zampini 501bd0d7efSStefano Zampini int main(int argc, char **argv) 511bd0d7efSStefano Zampini { 521bd0d7efSStefano Zampini DM dm; 531bd0d7efSStefano Zampini Mat M; 541bd0d7efSStefano Zampini Vec llM, lM, ones, work; 551bd0d7efSStefano Zampini PetscReal rtol = PETSC_SMALL, atol = 0.0; 561bd0d7efSStefano Zampini PetscScalar vals[6]; 571bd0d7efSStefano Zampini PetscInt dim; 581bd0d7efSStefano Zampini PetscBool amr = PETSC_FALSE; 591bd0d7efSStefano Zampini 601bd0d7efSStefano Zampini PetscFunctionBeginUser; 611bd0d7efSStefano Zampini PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 621bd0d7efSStefano Zampini PetscCall(PetscOptionsGetBool(NULL, NULL, "-amr", &amr, NULL)); 631bd0d7efSStefano Zampini PetscCall(DMCreate(PETSC_COMM_WORLD, &dm)); 641bd0d7efSStefano Zampini PetscCall(DMSetType(dm, DMPLEX)); 651bd0d7efSStefano Zampini PetscCall(DMSetFromOptions(dm)); 661bd0d7efSStefano Zampini PetscCall(DMGetDimension(dm, &dim)); 671bd0d7efSStefano Zampini if (amr) { 681bd0d7efSStefano Zampini DM dmConv; 691bd0d7efSStefano Zampini 701bd0d7efSStefano Zampini PetscCall(DMConvert(dm, dim == 2 ? DMP4EST : DMP8EST, &dmConv)); 711bd0d7efSStefano Zampini PetscCall(DMDestroy(&dm)); 721bd0d7efSStefano Zampini dm = dmConv; 731bd0d7efSStefano Zampini PetscCall(DMSetFromOptions(dm)); 741bd0d7efSStefano Zampini PetscCall(DMSetUp(dm)); 751bd0d7efSStefano Zampini } 761bd0d7efSStefano Zampini PetscCall(DMViewFromOptions(dm, NULL, "-mesh_view")); 771bd0d7efSStefano Zampini PetscCall(SetupDiscretization(dm)); 781bd0d7efSStefano Zampini PetscCall(DMGetGlobalVector(dm, &ones)); 791bd0d7efSStefano Zampini PetscCall(DMGetGlobalVector(dm, &work)); 801bd0d7efSStefano Zampini PetscCall(VecSet(ones, 1.0)); 811bd0d7efSStefano Zampini PetscCall(DMCreateMassMatrixLumped(dm, &llM, &lM)); 821bd0d7efSStefano Zampini PetscCall(VecViewFromOptions(llM, NULL, "-local_lumped_mass_view")); 831bd0d7efSStefano Zampini PetscCall(VecViewFromOptions(lM, NULL, "-lumped_mass_view")); 841bd0d7efSStefano Zampini PetscCall(DMCreateMassMatrix(dm, NULL, &M)); 851bd0d7efSStefano Zampini PetscCall(MatViewFromOptions(M, NULL, "-mass_view")); 861bd0d7efSStefano Zampini PetscCall(MatMult(M, ones, work)); 871bd0d7efSStefano Zampini PetscCall(VecViewFromOptions(work, NULL, "-mass_rowsum_view")); 881bd0d7efSStefano Zampini PetscCall(VecSum(work, &vals[3])); 891bd0d7efSStefano Zampini PetscCall(VecSet(work, 0.0)); 901bd0d7efSStefano Zampini PetscCall(DMLocalToGlobal(dm, llM, ADD_VALUES, work)); 911bd0d7efSStefano Zampini PetscCall(VecSum(work, &vals[4])); 921bd0d7efSStefano Zampini PetscCall(VecSum(lM, &vals[5])); 931bd0d7efSStefano Zampini PetscCall(DMPlexComputeIntegralFEM(dm, ones, vals, NULL)); 941bd0d7efSStefano Zampini CheckVals(vals[0], vals[1], rtol, atol, "Error volume"); 951bd0d7efSStefano Zampini CheckVals((3 + dim) * vals[0], vals[3], rtol, atol, "Error mass"); 961bd0d7efSStefano Zampini CheckVals((3 + dim) * vals[0], vals[4], rtol, atol, "Error local lumped mass"); 971bd0d7efSStefano Zampini CheckVals((3 + dim) * vals[0], vals[5], rtol, atol, "Error lumped mass"); 981bd0d7efSStefano Zampini PetscCall(MatDestroy(&M)); 991bd0d7efSStefano Zampini PetscCall(VecDestroy(&lM)); 1001bd0d7efSStefano Zampini PetscCall(VecDestroy(&llM)); 1011bd0d7efSStefano Zampini PetscCall(DMRestoreGlobalVector(dm, &work)); 1021bd0d7efSStefano Zampini PetscCall(DMRestoreGlobalVector(dm, &ones)); 1031bd0d7efSStefano Zampini PetscCall(DMDestroy(&dm)); 1041bd0d7efSStefano Zampini PetscCall(PetscFinalize()); 1051bd0d7efSStefano Zampini return 0; 1061bd0d7efSStefano Zampini } 1071bd0d7efSStefano Zampini 1081bd0d7efSStefano Zampini /*TEST 1091bd0d7efSStefano Zampini 1101bd0d7efSStefano Zampini test: 1111bd0d7efSStefano Zampini suffix: 0 1123886731fSPierre Jolivet output_file: output/empty.out 1131bd0d7efSStefano Zampini nsize: {{1 2}} 1141bd0d7efSStefano Zampini args: -dm_plex_box_faces 3,3,3 -dm_plex_dim {{2 3}} -dm_plex_simplex 0 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -petscpartitioner_type simple 1151bd0d7efSStefano Zampini 1161bd0d7efSStefano Zampini test: 1171bd0d7efSStefano Zampini TODO: broken 1181bd0d7efSStefano Zampini suffix: 0_amr 1193886731fSPierre Jolivet output_file: output/empty.out 1201bd0d7efSStefano Zampini requires: p4est 1211bd0d7efSStefano Zampini nsize: {{1 2}} 1221bd0d7efSStefano Zampini args: -dm_plex_dim {{2 3}} -dm_plex_simplex 0 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -amr -dm_p4est_refine_pattern hash -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 2 1231bd0d7efSStefano Zampini 1241bd0d7efSStefano Zampini test: 1251bd0d7efSStefano Zampini suffix: 1 1263886731fSPierre Jolivet output_file: output/empty.out 1271bd0d7efSStefano Zampini requires: triangle 1281bd0d7efSStefano Zampini nsize: {{1 2}} 1291bd0d7efSStefano Zampini args: -dm_plex_box_faces 3,3 -dm_plex_dim 2 -dm_plex_simplex 1 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -petscpartitioner_type simple 1301bd0d7efSStefano Zampini 1311bd0d7efSStefano Zampini test: 1321bd0d7efSStefano Zampini suffix: 2 1333886731fSPierre Jolivet output_file: output/empty.out 1341bd0d7efSStefano Zampini requires: ctetgen 1351bd0d7efSStefano Zampini nsize: {{1 2}} 1361bd0d7efSStefano Zampini args: -dm_plex_box_faces 3,3,3 -dm_plex_dim 3 -dm_plex_simplex 1 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -petscpartitioner_type simple 1371bd0d7efSStefano Zampini 1381bd0d7efSStefano Zampini TEST*/ 139