Lines Matching +full:- +full:- +full:prefix

1 // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
4 // SPDX-License-Identifier: BSD-2-Clause
112 (*ctx)->dm_x = dm_x; in OperatorApplyContextCreate()
114 (*ctx)->dm_y = dm_y; in OperatorApplyContextCreate()
117 (*ctx)->X_loc = X_loc; in OperatorApplyContextCreate()
119 (*ctx)->Y_loc = Y_loc; in OperatorApplyContextCreate()
122 if (x_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(x_ceed, &(*ctx)->x_ceed)); in OperatorApplyContextCreate()
123 else PetscCallCeed(ceed, CeedVectorCreate(ceed, x_size, &(*ctx)->x_ceed)); in OperatorApplyContextCreate()
125 if (y_ceed) PetscCallCeed(ceed, CeedVectorReferenceCopy(y_ceed, &(*ctx)->y_ceed)); in OperatorApplyContextCreate()
126 else PetscCallCeed(ceed, CeedVectorCreate(ceed, y_size, &(*ctx)->y_ceed)); in OperatorApplyContextCreate()
128 PetscCallCeed(ceed, CeedOperatorReferenceCopy(op_apply, &(*ctx)->op)); in OperatorApplyContextCreate()
129 PetscCallCeed(ceed, CeedReferenceCopy(ceed, &(*ctx)->ceed)); in OperatorApplyContextCreate()
141 Ceed ceed = ctx->ceed; in OperatorApplyContextDestroy()
144 PetscCall(DMDestroy(&ctx->dm_x)); in OperatorApplyContextDestroy()
145 PetscCall(DMDestroy(&ctx->dm_y)); in OperatorApplyContextDestroy()
146 PetscCall(VecDestroy(&ctx->X_loc)); in OperatorApplyContextDestroy()
147 PetscCall(VecDestroy(&ctx->Y_loc)); in OperatorApplyContextDestroy()
150 PetscCallCeed(ceed, CeedVectorDestroy(&ctx->x_ceed)); in OperatorApplyContextDestroy()
151 PetscCallCeed(ceed, CeedVectorDestroy(&ctx->y_ceed)); in OperatorApplyContextDestroy()
152 PetscCallCeed(ceed, CeedOperatorDestroy(&ctx->op)); in OperatorApplyContextDestroy()
153 PetscCallCeed(ceed, CeedDestroy(&ctx->ceed)); in OperatorApplyContextDestroy()
171 …* For example, if statitics are being store at quadrature points, a `DM`-created `Vec` will not ha…
219 Ceed ceed = ctx->ceed; in ApplyCeedOperator_Core()
222 if (X) PetscCall(DMGlobalToLocal(ctx->dm_x, X, INSERT_VALUES, X_loc)); in ApplyCeedOperator_Core()
229 …if (use_apply_add) PetscCallCeed(ceed, CeedOperatorApplyAdd(ctx->op, x_ceed, y_ceed, CEED_REQUEST_… in ApplyCeedOperator_Core()
230 else PetscCallCeed(ceed, CeedOperatorApply(ctx->op, x_ceed, y_ceed, CEED_REQUEST_IMMEDIATE)); in ApplyCeedOperator_Core()
234 if (X_loc) PetscCall(VecReadCeedToPetsc(ctx->x_ceed, x_mem_type, X_loc)); in ApplyCeedOperator_Core()
236 if (Y_loc) PetscCall(VecCeedToPetsc(ctx->y_ceed, y_mem_type, Y_loc)); in ApplyCeedOperator_Core()
237 if (Y) PetscCall(DMLocalToGlobal(ctx->dm_y, Y_loc, ADD_VALUES, Y)); in ApplyCeedOperator_Core()
242 Vec X_loc = ctx->X_loc, Y_loc = ctx->Y_loc; in ApplyCeedOperatorGlobalToGlobal()
248 if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); in ApplyCeedOperatorGlobalToGlobal()
249 if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); in ApplyCeedOperatorGlobalToGlobal()
251 PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false)); in ApplyCeedOperatorGlobalToGlobal()
254 if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); in ApplyCeedOperatorGlobalToGlobal()
255 if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); in ApplyCeedOperatorGlobalToGlobal()
260 Vec Y_loc = ctx->Y_loc; in ApplyCeedOperatorLocalToGlobal()
266 if (!ctx->Y_loc) PetscCall(DMGetLocalVector(ctx->dm_y, &Y_loc)); in ApplyCeedOperatorLocalToGlobal()
268 PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, Y, ctx, false)); in ApplyCeedOperatorLocalToGlobal()
271 if (!ctx->Y_loc) PetscCall(DMRestoreLocalVector(ctx->dm_y, &Y_loc)); in ApplyCeedOperatorLocalToGlobal()
276 Vec X_loc = ctx->X_loc; in ApplyCeedOperatorGlobalToLocal()
280 if (!ctx->X_loc) PetscCall(DMGetLocalVector(ctx->dm_x, &X_loc)); in ApplyCeedOperatorGlobalToLocal()
282 PetscCall(ApplyCeedOperator_Core(X, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false)); in ApplyCeedOperatorGlobalToLocal()
285 if (!ctx->X_loc) PetscCall(DMRestoreLocalVector(ctx->dm_x, &X_loc)); in ApplyCeedOperatorGlobalToLocal()
291 PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, false)); in ApplyCeedOperatorLocalToLocal()
297 PetscCall(ApplyCeedOperator_Core(NULL, X_loc, ctx->x_ceed, ctx->y_ceed, Y_loc, NULL, ctx, true)); in ApplyAddCeedOperatorLocalToLocal()
304 …* Uses command-line flag with `ksp`'s prefix to determine if mat_ceed should be used directly or w…
324 …PetscCall(PetscOptionsBool("-matceed_assemble_amat", "Assemble the A matrix for KSP solve", NULL, … in CreateSolveOperatorsFromMatCeed()
378 // -----------------------------------------------------------------------------