1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 /// @file 5 /// Utility functions for setting up Taylor-Green Vortex 6 7 #include "../qfunctions/taylorgreen.h" 8 9 #include <navierstokes.h> 10 11 PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx) { 12 Honee honee = *(Honee *)ctx; 13 Ceed ceed = honee->ceed; 14 MPI_Comm comm = honee->comm; 15 TaylorGreenContext taylorgreen_ctx; 16 CeedQFunctionContext taylorgreen_qfctx; 17 SetupContext setup_ctx; 18 PetscScalar background_velocity[3] = {0}; 19 PetscInt background_velocity_size = PETSC_STATIC_ARRAY_LENGTH(background_velocity); 20 21 PetscFunctionBeginUser; 22 PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx)); 23 24 PetscOptionsBegin(comm, NULL, "Options for TAYLOR-GREEN problem", NULL); 25 PetscCall(PetscOptionsScalarArray("-taylorgreen_background_velocity", "Background velocity to add to Taylor-Green initial condition", NULL, 26 background_velocity, &background_velocity_size, NULL)); 27 PetscOptionsEnd(); 28 29 PetscCall(PetscNew(&taylorgreen_ctx)); 30 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->ics.qfctx, CEED_MEM_HOST, &setup_ctx)); 31 taylorgreen_ctx->reference = setup_ctx->reference; 32 taylorgreen_ctx->gas = setup_ctx->gas; 33 taylorgreen_ctx->lx = setup_ctx->lx; 34 taylorgreen_ctx->ly = setup_ctx->ly; 35 taylorgreen_ctx->lz = setup_ctx->lz; 36 PetscCall(PetscArraycpy(taylorgreen_ctx->u, background_velocity, background_velocity_size)); 37 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->ics.qfctx, &setup_ctx)); 38 39 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &taylorgreen_qfctx)); 40 PetscCallCeed(ceed, CeedQFunctionContextSetData(taylorgreen_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*taylorgreen_ctx), taylorgreen_ctx)); 41 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(taylorgreen_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 42 43 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 44 problem->ics.qf_func_ptr = ICsTaylorGreen; 45 problem->ics.qf_loc = ICsTaylorGreen_loc; 46 problem->ics.qfctx = taylorgreen_qfctx; 47 PetscFunctionReturn(PETSC_SUCCESS); 48 } 49