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