xref: /libCEED/examples/fluids/problems/gaussianwave.c (revision ee4ca9cbfe2be39196684117442f3ce8d00b6506)
1530ad8c4SKenneth E. Jansen // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2530ad8c4SKenneth E. Jansen // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3530ad8c4SKenneth E. Jansen //
4530ad8c4SKenneth E. Jansen // SPDX-License-Identifier: BSD-2-Clause
5530ad8c4SKenneth E. Jansen //
6530ad8c4SKenneth E. Jansen // This file is part of CEED:  http://github.com/ceed
7530ad8c4SKenneth E. Jansen 
8530ad8c4SKenneth E. Jansen /// @file
9530ad8c4SKenneth E. Jansen /// Utility functions for setting up Gaussian Wave problem
10530ad8c4SKenneth E. Jansen 
11530ad8c4SKenneth E. Jansen #include "../qfunctions/gaussianwave.h"
12530ad8c4SKenneth E. Jansen 
1349aac155SJeremy L Thompson #include <ceed.h>
1449aac155SJeremy L Thompson #include <petscdm.h>
1549aac155SJeremy L Thompson 
16530ad8c4SKenneth E. Jansen #include "../navierstokes.h"
17530ad8c4SKenneth E. Jansen #include "../qfunctions/freestream_bc_type.h"
18530ad8c4SKenneth E. Jansen 
19530ad8c4SKenneth E. Jansen PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc) {
20530ad8c4SKenneth E. Jansen   User                     user = *(User *)ctx;
21530ad8c4SKenneth E. Jansen   MPI_Comm                 comm = PETSC_COMM_WORLD;
22530ad8c4SKenneth E. Jansen   GaussianWaveContext      gausswave_ctx;
23530ad8c4SKenneth E. Jansen   FreestreamContext        freestream_ctx;
24530ad8c4SKenneth E. Jansen   NewtonianIdealGasContext newtonian_ig_ctx;
25530ad8c4SKenneth E. Jansen   CeedQFunctionContext     gausswave_context;
26530ad8c4SKenneth E. Jansen 
27530ad8c4SKenneth E. Jansen   PetscFunctionBeginUser;
28530ad8c4SKenneth E. Jansen   PetscCall(NS_NEWTONIAN_IG(problem, dm, ctx, bc));
29530ad8c4SKenneth E. Jansen 
30530ad8c4SKenneth E. Jansen   switch (user->phys->state_var) {
31530ad8c4SKenneth E. Jansen     case STATEVAR_CONSERVATIVE:
32530ad8c4SKenneth E. Jansen       problem->ics.qfunction     = IC_GaussianWave_Conserv;
33530ad8c4SKenneth E. Jansen       problem->ics.qfunction_loc = IC_GaussianWave_Conserv_loc;
34530ad8c4SKenneth E. Jansen       break;
35530ad8c4SKenneth E. Jansen     case STATEVAR_PRIMITIVE:
36530ad8c4SKenneth E. Jansen       problem->ics.qfunction     = IC_GaussianWave_Prim;
37530ad8c4SKenneth E. Jansen       problem->ics.qfunction_loc = IC_GaussianWave_Prim_loc;
38530ad8c4SKenneth E. Jansen       break;
39530ad8c4SKenneth E. Jansen   }
40530ad8c4SKenneth E. Jansen 
41530ad8c4SKenneth E. Jansen   // -- Option Defaults
42530ad8c4SKenneth E. Jansen   CeedScalar epicenter[3] = {0.};   // m
43530ad8c4SKenneth E. Jansen   CeedScalar width        = 0.002;  // m
44530ad8c4SKenneth E. Jansen   CeedScalar amplitude    = 0.1;    // -
45530ad8c4SKenneth E. Jansen 
46530ad8c4SKenneth E. Jansen   PetscOptionsBegin(comm, NULL, "Options for GAUSSIAN_WAVE problem", NULL);
47530ad8c4SKenneth E. Jansen   PetscInt narray = 3;
48530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalarArray("-epicenter", "Coordinates of center of wave", NULL, epicenter, &narray, NULL));
49530ad8c4SKenneth E. Jansen   PetscCheck(narray == 3, comm, PETSC_ERR_ARG_SIZ, "-epicenter should recieve array of size 3, instead recieved size %" PetscInt_FMT ".", narray);
50530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-width", "Width parameter for perturbation size", NULL, width, &width, NULL));
51530ad8c4SKenneth E. Jansen   PetscCall(PetscOptionsScalar("-amplitude", "Amplitude of the perturbation", NULL, amplitude, &amplitude, NULL));
52530ad8c4SKenneth E. Jansen   PetscOptionsEnd();
53530ad8c4SKenneth E. Jansen 
54530ad8c4SKenneth E. Jansen   width *= user->units->meter;
55530ad8c4SKenneth E. Jansen   for (int i = 0; i < 3; i++) epicenter[i] *= user->units->meter;
56530ad8c4SKenneth E. Jansen 
57530ad8c4SKenneth E. Jansen   // -- Set gausswave_ctx struct values
58530ad8c4SKenneth E. Jansen   PetscCall(PetscCalloc1(1, &gausswave_ctx));
59530ad8c4SKenneth E. Jansen   CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
60530ad8c4SKenneth E. Jansen   CeedQFunctionContextGetData(problem->apply_freestream.qfunction_context, CEED_MEM_HOST, &freestream_ctx);
61530ad8c4SKenneth E. Jansen 
62530ad8c4SKenneth E. Jansen   gausswave_ctx->amplitude = amplitude;
63530ad8c4SKenneth E. Jansen   gausswave_ctx->width     = width;
64530ad8c4SKenneth E. Jansen   gausswave_ctx->S_infty   = freestream_ctx->S_infty;
65530ad8c4SKenneth E. Jansen   gausswave_ctx->newt_ctx  = *newtonian_ig_ctx;
66530ad8c4SKenneth E. Jansen   PetscCall(PetscArraycpy(gausswave_ctx->epicenter, epicenter, 3));
67530ad8c4SKenneth E. Jansen 
68530ad8c4SKenneth E. Jansen   CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
69530ad8c4SKenneth E. Jansen   CeedQFunctionContextRestoreData(problem->apply_freestream.qfunction_context, &freestream_ctx);
70530ad8c4SKenneth E. Jansen 
71530ad8c4SKenneth E. Jansen   CeedQFunctionContextCreate(user->ceed, &gausswave_context);
72530ad8c4SKenneth E. Jansen   CeedQFunctionContextSetData(gausswave_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*gausswave_ctx), gausswave_ctx);
73530ad8c4SKenneth E. Jansen   CeedQFunctionContextSetDataDestroy(gausswave_context, CEED_MEM_HOST, FreeContextPetsc);
74530ad8c4SKenneth E. Jansen   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
75530ad8c4SKenneth E. Jansen   problem->ics.qfunction_context = gausswave_context;
76530ad8c4SKenneth E. Jansen 
77*ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
78530ad8c4SKenneth E. Jansen }
79