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 Gaussian Wave problem 6 7 #include "../qfunctions/gaussianwave.h" 8 9 #include <ceed.h> 10 #include <petscdm.h> 11 12 #include <navierstokes.h> 13 #include "../qfunctions/bc_freestream_type.h" 14 15 PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc) { 16 Honee honee = *(Honee *)ctx; 17 MPI_Comm comm = honee->comm; 18 Ceed ceed = honee->ceed; 19 GaussianWaveContext gausswave_ctx; 20 FreestreamContext freestream_ctx; 21 NewtonianIdealGasContext newtonian_ig_ctx; 22 CeedQFunctionContext gausswave_qfctx; 23 24 PetscFunctionBeginUser; 25 PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc)); 26 27 switch (honee->phys->state_var) { 28 case STATEVAR_CONSERVATIVE: 29 problem->ics.qf_func_ptr = IC_GaussianWave_Conserv; 30 problem->ics.qf_loc = IC_GaussianWave_Conserv_loc; 31 break; 32 case STATEVAR_PRIMITIVE: 33 problem->ics.qf_func_ptr = IC_GaussianWave_Prim; 34 problem->ics.qf_loc = IC_GaussianWave_Prim_loc; 35 break; 36 case STATEVAR_ENTROPY: 37 problem->ics.qf_func_ptr = IC_GaussianWave_Entropy; 38 problem->ics.qf_loc = IC_GaussianWave_Entropy_loc; 39 break; 40 } 41 42 // -- Option Defaults 43 CeedScalar epicenter[3] = {0.}; // m 44 CeedScalar width = 0.002; // m 45 CeedScalar amplitude = 0.1; // - 46 47 PetscOptionsBegin(comm, NULL, "Options for GAUSSIAN_WAVE problem", NULL); 48 PetscInt narray = 3; 49 PetscCall(PetscOptionsScalarArray("-epicenter", "Coordinates of center of wave", NULL, epicenter, &narray, NULL)); 50 PetscCheck(narray == 3, comm, PETSC_ERR_ARG_SIZ, "-epicenter should recieve array of size 3, instead recieved size %" PetscInt_FMT ".", narray); 51 PetscCall(PetscOptionsScalar("-width", "Width parameter for perturbation size", NULL, width, &width, NULL)); 52 PetscCall(PetscOptionsScalar("-amplitude", "Amplitude of the perturbation", NULL, amplitude, &litude, NULL)); 53 PetscOptionsEnd(); 54 55 width *= honee->units->meter; 56 for (int i = 0; i < 3; i++) epicenter[i] *= honee->units->meter; 57 58 // -- Set gausswave_ctx struct values 59 PetscCall(PetscCalloc1(1, &gausswave_ctx)); 60 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 61 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_freestream.qfctx, CEED_MEM_HOST, &freestream_ctx)); 62 63 gausswave_ctx->amplitude = amplitude; 64 gausswave_ctx->width = width; 65 gausswave_ctx->S_infty = freestream_ctx->S_infty; 66 gausswave_ctx->newt_ctx = *newtonian_ig_ctx; 67 PetscCall(PetscArraycpy(gausswave_ctx->epicenter, epicenter, 3)); 68 69 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 70 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_freestream.qfctx, &freestream_ctx)); 71 72 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &gausswave_qfctx)); 73 PetscCallCeed(ceed, CeedQFunctionContextSetData(gausswave_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx)); 74 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(gausswave_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 75 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfctx)); 76 problem->ics.qfctx = gausswave_qfctx; 77 PetscFunctionReturn(PETSC_SUCCESS); 78 } 79