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