1 // Copyright (c) 2017-2022, 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/freestream_bc_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 } 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 *= user->units->meter; 56 for (int i = 0; i < 3; i++) epicenter[i] *= user->units->meter; 57 58 // -- Set gausswave_ctx struct values 59 PetscCall(PetscCalloc1(1, &gausswave_ctx)); 60 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 61 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_freestream.qfunction_context, 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.qfunction_context, &newtonian_ig_ctx)); 70 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_freestream.qfunction_context, &freestream_ctx)); 71 72 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &gausswave_context)); 73 PetscCallCeed(ceed, CeedQFunctionContextSetData(gausswave_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx)); 74 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(gausswave_context, CEED_MEM_HOST, FreeContextPetsc)); 75 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&problem->ics.qfunction_context)); 76 problem->ics.qfunction_context = gausswave_context; 77 PetscFunctionReturn(PETSC_SUCCESS); 78 } 79