1c4762a1bSJed Brown static char help[] = "Tests periodic boundary conditions for DMDA1d with periodic boundary conditions.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmda.h> 4c4762a1bSJed Brown 5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 6d71ae5a4SJacob Faibussowitsch { 7c4762a1bSJed Brown DM da; 8c4762a1bSJed Brown Mat A; 9c4762a1bSJed Brown const PetscInt dfill[4] = {0, 1, 0, 1}, ofill[4] = {0, 1, 1, 0}; 10c4762a1bSJed Brown 11327415f7SBarry Smith PetscFunctionBeginUser; 12*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 139566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, 7, 2, 1, NULL, &da)); 149566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 159566063dSJacob Faibussowitsch PetscCall(DMDASetBlockFills(da, dfill, ofill)); 169566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 179566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &A)); 189566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 199566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_INFO)); 209566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD)); 219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 229566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 239566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 24b122ec5aSJacob Faibussowitsch return 0; 25c4762a1bSJed Brown } 26c4762a1bSJed Brown 27c4762a1bSJed Brown /*TEST 28c4762a1bSJed Brown 29c4762a1bSJed Brown test: 30c4762a1bSJed Brown 31c4762a1bSJed Brown test: 32c4762a1bSJed Brown suffix: 2 33c4762a1bSJed Brown nsize: 2 34c4762a1bSJed Brown 35c4762a1bSJed Brown TEST*/ 36