Lines Matching +full:- +full:f

1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
19 …comps Number of components that makes up the diffusive flux (e.g. 1 for scalar advection-diffusion)
24 … height = 0, dm_field = 0, dim, degree = honee->app_ctx->degree, q_extra = honee->app_ctx- in DivDiffFluxProjectionCreate()
41 PetscCall(DMClone(honee->dm, &projection->dm)); in DivDiffFluxProjectionCreate()
42 PetscCall(DMSetMatrixPreallocateSkip(projection->dm, PETSC_TRUE)); in DivDiffFluxProjectionCreate()
43 PetscCall(DMGetDimension(projection->dm, &dim)); in DivDiffFluxProjectionCreate()
44 switch (diff_flux_proj_->method) { in DivDiffFluxProjectionCreate()
46 projection->num_comp = diff_flux_proj_->num_diff_flux_comps; in DivDiffFluxProjectionCreate()
47 PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DivDiffFluxProj")); in DivDiffFluxProjectionCreate()
48 …ByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); in DivDiffFluxProjectionCreate()
50 …PetscCall(DMPlexCeedElemRestrictionCreate(honee->ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DE… in DivDiffFluxProjectionCreate()
51 &diff_flux_proj_->elem_restr_div_diff_flux)); in DivDiffFluxProjectionCreate()
52 PetscCallCeed(honee->ceed, in DivDiffFluxProjectionCreate()
53 …CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_d… in DivDiffFluxProjectionCreate()
54 …PetscCall(DMPlexCeedBasisCreate(honee->ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALU… in DivDiffFluxProjectionCreate()
55 &diff_flux_proj_->basis_div_diff_flux)); in DivDiffFluxProjectionCreate()
56 diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_INTERP; in DivDiffFluxProjectionCreate()
58 { // Create face labels on projection->dm for boundary integrals in DivDiffFluxProjectionCreate()
62 PetscCall(DMGetLabel(honee->dm, "Face Sets", &face_sets_label)); in DivDiffFluxProjectionCreate()
63 …PetscCall(DMLabelCreateGlobalValueArray(honee->dm, face_sets_label, &num_face_set_values, &face_se… in DivDiffFluxProjectionCreate()
64 for (PetscInt f = 0; f < num_face_set_values; f++) { in DivDiffFluxProjectionCreate() local
68 … PetscCall(DMPlexCreateFaceLabel(honee->dm, face_set_values[f], &face_orientation_label_name)); in DivDiffFluxProjectionCreate()
69 PetscCall(DMGetLabel(honee->dm, face_orientation_label_name, &face_orientation_label)); in DivDiffFluxProjectionCreate()
70 PetscCall(DMAddLabel(projection->dm, face_orientation_label)); in DivDiffFluxProjectionCreate()
77 projection->num_comp = diff_flux_proj_->num_diff_flux_comps * dim; in DivDiffFluxProjectionCreate()
78 PetscCall(PetscObjectSetName((PetscObject)projection->dm, "DiffFluxProj")); in DivDiffFluxProjectionCreate()
79 …ByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &projection->num_comp, projection->dm)); in DivDiffFluxProjectionCreate()
81 …PetscCall(DMPlexCeedElemRestrictionQDataCreate(honee->ceed, projection->dm, DMLABEL_DEFAULT, DMLAB… in DivDiffFluxProjectionCreate()
82 … diff_flux_proj_->num_diff_flux_comps, &diff_flux_proj_->elem_restr_div_diff_flux)); in DivDiffFluxProjectionCreate()
83 PetscCallCeed(honee->ceed, in DivDiffFluxProjectionCreate()
84 …CeedElemRestrictionCreateVector(diff_flux_proj_->elem_restr_div_diff_flux, &diff_flux_proj_->div_d… in DivDiffFluxProjectionCreate()
85 diff_flux_proj_->basis_div_diff_flux = CEED_BASIS_NONE; in DivDiffFluxProjectionCreate()
86 diff_flux_proj_->eval_mode_div_diff_flux = CEED_EVAL_NONE; in DivDiffFluxProjectionCreate()
89 …SETERRQ(honee->comm, PETSC_ERR_ARG_WRONG, "Should not reach here with div_diff_flux_projection_met… in DivDiffFluxProjectionCreate()
108 Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); in DivDiffFluxProjectionGetOperatorFieldData()
111 …if (elem_restr) PetscCallCeed(ceed, CeedElemRestrictionReferenceCopy(diff_flux_proj->elem_restr_di… in DivDiffFluxProjectionGetOperatorFieldData()
112 …if (basis) PetscCallCeed(ceed, CeedBasisReferenceCopy(diff_flux_proj->basis_div_diff_flux, basis)); in DivDiffFluxProjectionGetOperatorFieldData()
113 …if (vector) PetscCallCeed(ceed, CeedVectorReferenceCopy(diff_flux_proj->div_diff_flux_ceed, vector… in DivDiffFluxProjectionGetOperatorFieldData()
114 if (eval_mode) *eval_mode = diff_flux_proj->eval_mode_div_diff_flux; in DivDiffFluxProjectionGetOperatorFieldData()
125 Ceed ceed = honee->ceed; in DivDiffFluxProjectionSetup_Direct()
126 NodalProjectionData projection = diff_flux_proj->projection; in DivDiffFluxProjectionSetup_Direct()
127 MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); in DivDiffFluxProjectionSetup_Direct()
133 PetscCheck(diff_flux_proj->CreateRHSOperator_Direct, comm, PETSC_ERR_ARG_WRONGSTATE, in DivDiffFluxProjectionSetup_Direct()
135 PetscCall(diff_flux_proj->CreateRHSOperator_Direct(honee, diff_flux_proj, &op_rhs)); in DivDiffFluxProjectionSetup_Direct()
136 PetscCall(DMCreateLocalVector(projection->dm, &diff_flux_proj->DivDiffFlux_loc)); in DivDiffFluxProjectionSetup_Direct()
137 diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; in DivDiffFluxProjectionSetup_Direct()
138 …PetscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, di… in DivDiffFluxProjectionSetup_Direct()
139 &projection->l2_rhs_ctx)); in DivDiffFluxProjectionSetup_Direct()
143 { // -- Build Mass operator in DivDiffFluxProjectionSetup_Direct()
152 …PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &… in DivDiffFluxProjectionSetup_Direct()
154 PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); in DivDiffFluxProjectionSetup_Direct()
160 { // -- Setup KSP for L^2 projection in DivDiffFluxProjectionSetup_Direct()
163 PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); in DivDiffFluxProjectionSetup_Direct()
165 PetscCall(KSPCreate(comm, &projection->ksp)); in DivDiffFluxProjectionSetup_Direct()
166 PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); in DivDiffFluxProjectionSetup_Direct()
169 PetscCall(KSPGetPC(projection->ksp, &pc)); in DivDiffFluxProjectionSetup_Direct()
172 PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); in DivDiffFluxProjectionSetup_Direct()
174 PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); in DivDiffFluxProjectionSetup_Direct()
195 Ceed ceed = honee->ceed; in DivDiffFluxProjectionSetup_Indirect()
196 NodalProjectionData projection = diff_flux_proj->projection; in DivDiffFluxProjectionSetup_Indirect()
201 MPI_Comm comm = PetscObjectComm((PetscObject)projection->dm); in DivDiffFluxProjectionSetup_Indirect()
207 …PetscCall(DMPlexCeedElemRestrictionCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_V… in DivDiffFluxProjectionSetup_Indirect()
208 …PetscCall(DMPlexCeedBasisCreate(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, heig… in DivDiffFluxProjectionSetup_Indirect()
209 …PetscCall(QDataGet(ceed, projection->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, &elem_restr_qd, &… in DivDiffFluxProjectionSetup_Indirect()
215 PetscCheck(diff_flux_proj->CreateRHSOperator_Indirect, comm, PETSC_ERR_ARG_WRONGSTATE, in DivDiffFluxProjectionSetup_Indirect()
217 PetscCall(diff_flux_proj->CreateRHSOperator_Indirect(honee, diff_flux_proj, &op_rhs)); in DivDiffFluxProjectionSetup_Indirect()
218 …tscCall(OperatorApplyContextCreate(honee->dm, projection->dm, ceed, op_rhs, NULL, NULL, NULL, NULL… in DivDiffFluxProjectionSetup_Indirect()
222 { // -- Build Mass operator in DivDiffFluxProjectionSetup_Indirect()
226 PetscCall(HoneeMassQFunctionCreate(ceed, projection->num_comp, q_data_size, &qf_mass)); in DivDiffFluxProjectionSetup_Indirect()
232 { // -- Setup KSP for L^2 projection in DivDiffFluxProjectionSetup_Indirect()
235 PetscCall(MatCreateCeed(projection->dm, projection->dm, op_mass, NULL, &mat_mass)); in DivDiffFluxProjectionSetup_Indirect()
237 PetscCall(KSPCreate(comm, &projection->ksp)); in DivDiffFluxProjectionSetup_Indirect()
238 PetscCall(KSPSetOptionsPrefix(projection->ksp, "div_diff_flux_projection_")); in DivDiffFluxProjectionSetup_Indirect()
241 PetscCall(KSPGetPC(projection->ksp, &pc)); in DivDiffFluxProjectionSetup_Indirect()
244 PetscCall(KSPSetType(projection->ksp, KSPPREONLY)); in DivDiffFluxProjectionSetup_Indirect()
246 PetscCall(KSPSetFromOptions_WithMatCeed(projection->ksp, mat_mass)); in DivDiffFluxProjectionSetup_Indirect()
259 PetscCall(DMGetDimension(projection->dm, &dim)); in DivDiffFluxProjectionSetup_Indirect()
264 switch (diff_flux_proj->num_diff_flux_comps) { in DivDiffFluxProjectionSetup_Indirect()
272 switch (diff_flux_proj->num_diff_flux_comps) { in DivDiffFluxProjectionSetup_Indirect()
288 dim, diff_flux_proj->num_diff_flux_comps); in DivDiffFluxProjectionSetup_Indirect()
290 …PetscCallCeed(ceed, CeedQFunctionAddInput(qf_calc_divergence, "Grad F_diff", projection->num_comp … in DivDiffFluxProjectionSetup_Indirect()
292 …PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_calc_divergence, "Div F_diff", diff_flux_proj->num_d… in DivDiffFluxProjectionSetup_Indirect()
298 diff_flux_proj->div_diff_flux_ceed)); in DivDiffFluxProjectionSetup_Indirect()
300 …PetscCall(OperatorApplyContextCreate(projection->dm, NULL, ceed, op_calc_divergence, NULL, CEED_VE… in DivDiffFluxProjectionSetup_Indirect()
301 &diff_flux_proj->calc_div_diff_flux)); in DivDiffFluxProjectionSetup_Indirect()
322 switch (honee->app_ctx->divFdiffproj_method) { in DivDiffFluxProjectionSetup()
330 …SETERRQ(PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_ARG_WRONG, "Should not reach here with … in DivDiffFluxProjectionSetup()
331 DivDiffFluxProjectionMethods[honee->app_ctx->divFdiffproj_method]); in DivDiffFluxProjectionSetup()
346 NodalProjectionData projection = diff_flux_proj->projection; in DivDiffFluxProjectionApply()
350 switch (diff_flux_proj->method) { in DivDiffFluxProjectionApply()
354 PetscCall(DMGetGlobalVector(projection->dm, &DivDiffFlux)); in DivDiffFluxProjectionApply()
355 PetscCall(DMGetGlobalVector(projection->dm, &RHS)); in DivDiffFluxProjectionApply()
356 if (diff_flux_proj->ceed_vec_has_array) { in DivDiffFluxProjectionApply()
357 …cCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, … in DivDiffFluxProjectionApply()
358 diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; in DivDiffFluxProjectionApply()
360 PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); in DivDiffFluxProjectionApply()
361 PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); in DivDiffFluxProjectionApply()
364 // Run PCApply manually if using ksp_type preonly -pc_type jacobi in DivDiffFluxProjectionApply()
366 … // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix in DivDiffFluxProjectionApply()
369 PetscCall(KSPGetPC(projection->ksp, &pc)); in DivDiffFluxProjectionApply()
370 PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); in DivDiffFluxProjectionApply()
373 else PetscCall(KSPSolve(projection->ksp, RHS, DivDiffFlux)); in DivDiffFluxProjectionApply()
375 PetscCall(VecViewFromOptions(DivDiffFlux, NULL, "-div_diff_flux_projection_view")); in DivDiffFluxProjectionApply()
377 …PetscCall(DMGlobalToLocal(projection->dm, DivDiffFlux, INSERT_VALUES, diff_flux_proj->DivDiffFlux_… in DivDiffFluxProjectionApply()
378 …scCall(VecReadPetscToCeed(diff_flux_proj->DivDiffFlux_loc, &diff_flux_proj->DivDiffFlux_memtype, d… in DivDiffFluxProjectionApply()
379 diff_flux_proj->ceed_vec_has_array = PETSC_TRUE; in DivDiffFluxProjectionApply()
381 PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); in DivDiffFluxProjectionApply()
382 PetscCall(DMRestoreGlobalVector(projection->dm, &DivDiffFlux)); in DivDiffFluxProjectionApply()
388 PetscCall(DMGetGlobalVector(projection->dm, &DiffFlux)); in DivDiffFluxProjectionApply()
389 PetscCall(DMGetGlobalVector(projection->dm, &RHS)); in DivDiffFluxProjectionApply()
390 PetscCall(ApplyCeedOperatorLocalToGlobal(Q_loc, RHS, projection->l2_rhs_ctx)); in DivDiffFluxProjectionApply()
391 PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_rhs_view")); in DivDiffFluxProjectionApply()
394 // Run PCApply manually if using -ksp_type preonly -pc_type jacobi in DivDiffFluxProjectionApply()
396 … // See https://gitlab.com/petsc/petsc/-/merge_requests/8048 for more details and a possible fix in DivDiffFluxProjectionApply()
399 PetscCall(KSPGetPC(projection->ksp, &pc)); in DivDiffFluxProjectionApply()
400 PetscCall(PetscObjectTypeCompare((PetscObject)projection->ksp, KSPPREONLY, &ispreonly)); in DivDiffFluxProjectionApply()
403 else PetscCall(KSPSolve(projection->ksp, RHS, DiffFlux)); in DivDiffFluxProjectionApply()
405 PetscCall(VecViewFromOptions(DiffFlux, NULL, "-div_diff_flux_projection_view")); in DivDiffFluxProjectionApply()
407 PetscCall(ApplyCeedOperatorGlobalToLocal(DiffFlux, NULL, diff_flux_proj->calc_div_diff_flux)); in DivDiffFluxProjectionApply()
408 PetscCall(DMRestoreGlobalVector(projection->dm, &RHS)); in DivDiffFluxProjectionApply()
409 PetscCall(DMRestoreGlobalVector(projection->dm, &DiffFlux)); in DivDiffFluxProjectionApply()
412 …SETERRQ(PetscObjectComm((PetscObject)projection->dm), PETSC_ERR_ARG_WRONG, "Should not reach here … in DivDiffFluxProjectionApply()
413 DivDiffFluxProjectionMethods[diff_flux_proj->method]); in DivDiffFluxProjectionApply()
428 Ceed ceed = CeedVectorReturnCeed(diff_flux_proj->div_diff_flux_ceed); in DivDiffFluxProjectionDataDestroy()
430 PetscCall(NodalProjectionDataDestroy(&diff_flux_proj->projection)); in DivDiffFluxProjectionDataDestroy()
431 PetscCall(OperatorApplyContextDestroy(diff_flux_proj->calc_div_diff_flux)); in DivDiffFluxProjectionDataDestroy()
432 if (diff_flux_proj->ceed_vec_has_array) { in DivDiffFluxProjectionDataDestroy()
433 …cCall(VecReadCeedToPetsc(diff_flux_proj->div_diff_flux_ceed, diff_flux_proj->DivDiffFlux_memtype, … in DivDiffFluxProjectionDataDestroy()
434 diff_flux_proj->ceed_vec_has_array = PETSC_FALSE; in DivDiffFluxProjectionDataDestroy()
436 PetscCallCeed(ceed, CeedVectorDestroy(&diff_flux_proj->div_diff_flux_ceed)); in DivDiffFluxProjectionDataDestroy()
437 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&diff_flux_proj->elem_restr_div_diff_flux)); in DivDiffFluxProjectionDataDestroy()
438 PetscCallCeed(ceed, CeedBasisDestroy(&diff_flux_proj->basis_div_diff_flux)); in DivDiffFluxProjectionDataDestroy()
439 PetscCall(VecDestroy(&diff_flux_proj->DivDiffFlux_loc)); in DivDiffFluxProjectionDataDestroy()