135313810SToby Isaac const char help[] = "Test boundary condition insertion"; 235313810SToby Isaac 335313810SToby Isaac #include <petscdmplex.h> 435313810SToby Isaac 5*2a8381b2SBarry Smith static PetscErrorCode set_one(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], PetscCtx ctx) 6d71ae5a4SJacob Faibussowitsch { 735313810SToby Isaac PetscFunctionBegin; 835313810SToby Isaac bcval[0] = 1.; 93ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1035313810SToby Isaac } 1135313810SToby Isaac 12*2a8381b2SBarry Smith static PetscErrorCode set_two(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar bcval[], PetscCtx ctx) 13d71ae5a4SJacob Faibussowitsch { 1435313810SToby Isaac PetscFunctionBegin; 1535313810SToby Isaac bcval[0] = 2.; 163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1735313810SToby Isaac } 1835313810SToby Isaac 19d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 20d71ae5a4SJacob Faibussowitsch { 2135313810SToby Isaac DM dm; 2235313810SToby Isaac DMLabel label; 2335313810SToby Isaac PetscInt in_value = 1; 2435313810SToby Isaac PetscInt out_value = 3; 2535313810SToby Isaac PetscInt comps[] = {0}; 2635313810SToby Isaac PetscFE fe; 2735313810SToby Isaac Vec localVec; 2835313810SToby Isaac 29327415f7SBarry Smith PetscFunctionBeginUser; 3035313810SToby Isaac PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 3142108689Sksagiyam PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, 2, PETSC_FALSE, NULL, NULL, NULL, NULL, PETSC_TRUE, 0, PETSC_TRUE, &dm)); 3235313810SToby Isaac PetscCall(DMGetLabel(dm, "Face Sets", &label)); 3335313810SToby Isaac PetscCall(PetscFECreateLagrange(PETSC_COMM_WORLD, 2, 1, PETSC_FALSE, 1, PETSC_DETERMINE, &fe)); 3435313810SToby Isaac PetscCall(DMAddField(dm, NULL, (PetscObject)fe)); 3535313810SToby Isaac PetscCall(PetscFEDestroy(&fe)); 3635313810SToby Isaac PetscCall(DMCreateDS(dm)); 3757d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "inflow condition", label, 1, &in_value, 0, 1, comps, (PetscVoidFn *)set_one, NULL, NULL, NULL)); 3857d50842SBarry Smith PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "outflow condition", label, 1, &out_value, 0, 1, comps, (PetscVoidFn *)set_two, NULL, NULL, NULL)); 3935313810SToby Isaac PetscCall(DMCreateLocalVector(dm, &localVec)); 4035313810SToby Isaac PetscCall(VecSet(localVec, 0.)); 4135313810SToby Isaac PetscCall(DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localVec, 0.0, NULL, NULL, NULL)); 4235313810SToby Isaac PetscCall(VecView(localVec, NULL)); 4335313810SToby Isaac PetscCall(VecDestroy(&localVec)); 4435313810SToby Isaac PetscCall(DMDestroy(&dm)); 4535313810SToby Isaac PetscCall(PetscFinalize()); 4635313810SToby Isaac return 0; 4735313810SToby Isaac } 4835313810SToby Isaac 4935313810SToby Isaac /*TEST 5035313810SToby Isaac 5135313810SToby Isaac test: 5235313810SToby Isaac suffix: 0 5335313810SToby Isaac 5435313810SToby Isaac TEST*/ 55