xref: /libCEED/examples/fluids/problems/densitycurrent.c (revision 3d8e882215d238700cdceb37404f76ca7fa24eaa)
1*3d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*3d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
4*3d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
6*3d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Utility functions for setting up DENSITY_CURRENT
1077841947SLeila Ghaffari 
1177841947SLeila Ghaffari #include "../navierstokes.h"
1277841947SLeila Ghaffari #include "../qfunctions/densitycurrent.h"
1377841947SLeila Ghaffari 
141864f1c2SLeila Ghaffari PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *setup_ctx,
1577841947SLeila Ghaffari                                   void *ctx) {
1688b783a1SJames Wright 
1788b783a1SJames Wright   PetscInt ierr;
1888b783a1SJames Wright   ierr = NS_NEWTONIAN_IG(problem, dm, setup_ctx, ctx); CHKERRQ(ierr);
1977841947SLeila Ghaffari   SetupContext setup_context = *(SetupContext *)setup_ctx;
2077841947SLeila Ghaffari   User         user          = *(User *)ctx;
2177841947SLeila Ghaffari   MPI_Comm     comm          = PETSC_COMM_WORLD;
2277841947SLeila Ghaffari   PetscFunctionBeginUser;
2377841947SLeila Ghaffari 
2477841947SLeila Ghaffari   // ------------------------------------------------------
2577841947SLeila Ghaffari   //               SET UP DENSITY_CURRENT
2677841947SLeila Ghaffari   // ------------------------------------------------------
2777841947SLeila Ghaffari   problem->ics     = ICsDC;
2877841947SLeila Ghaffari   problem->ics_loc = ICsDC_loc;
2977841947SLeila Ghaffari   problem->bc      = Exact_DC;
3077841947SLeila Ghaffari 
3177841947SLeila Ghaffari   // ------------------------------------------------------
3277841947SLeila Ghaffari   //             Create the libCEED context
3377841947SLeila Ghaffari   // ------------------------------------------------------
3477841947SLeila Ghaffari   CeedScalar rc     = 1000.;   // m (Radius of bubble)
3577841947SLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0};
361864f1c2SLeila Ghaffari   PetscReal domain_min[3], domain_max[3], domain_size[3];
371864f1c2SLeila Ghaffari   ierr = DMGetBoundingBox(dm, domain_min, domain_max); CHKERRQ(ierr);
381864f1c2SLeila Ghaffari   for (int i=0; i<3; i++) domain_size[i] = domain_max[i] - domain_min[i];
3977841947SLeila Ghaffari 
4077841947SLeila Ghaffari   // ------------------------------------------------------
4177841947SLeila Ghaffari   //              Command line Options
4277841947SLeila Ghaffari   // ------------------------------------------------------
4377841947SLeila Ghaffari   ierr = PetscOptionsBegin(comm, NULL, "Options for DENSITY_CURRENT problem",
4477841947SLeila Ghaffari                            NULL); CHKERRQ(ierr);
4577841947SLeila Ghaffari   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
4677841947SLeila Ghaffari                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
471864f1c2SLeila Ghaffari   for (int i=0; i<3; i++) center[i] = .5*domain_size[i];
4877841947SLeila Ghaffari   PetscInt n = problem->dim;
4977841947SLeila Ghaffari   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
5077841947SLeila Ghaffari                                NULL, center, &n, NULL); CHKERRQ(ierr);
5177841947SLeila Ghaffari   n = problem->dim;
5277841947SLeila Ghaffari   ierr = PetscOptionsRealArray("-dc_axis",
5377841947SLeila Ghaffari                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
5477841947SLeila Ghaffari                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
5577841947SLeila Ghaffari   {
5677841947SLeila Ghaffari     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + PetscSqr(dc_axis[1]) +
5777841947SLeila Ghaffari                                    PetscSqr(dc_axis[2]));
5877841947SLeila Ghaffari     if (norm > 0) {
5977841947SLeila Ghaffari       for (int i=0; i<3; i++)  dc_axis[i] /= norm;
6077841947SLeila Ghaffari     }
6177841947SLeila Ghaffari   }
6277841947SLeila Ghaffari 
6377841947SLeila Ghaffari   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
6477841947SLeila Ghaffari 
6588b783a1SJames Wright   PetscScalar meter = user->units->meter;
6677841947SLeila Ghaffari   rc = fabs(rc) * meter;
6777841947SLeila Ghaffari   for (int i=0; i<3; i++) center[i] *= meter;
6877841947SLeila Ghaffari 
6977841947SLeila Ghaffari   setup_context->rc         = rc;
7077841947SLeila Ghaffari   setup_context->center[0]  = center[0];
7177841947SLeila Ghaffari   setup_context->center[1]  = center[1];
7277841947SLeila Ghaffari   setup_context->center[2]  = center[2];
7377841947SLeila Ghaffari   setup_context->dc_axis[0] = dc_axis[0];
7477841947SLeila Ghaffari   setup_context->dc_axis[1] = dc_axis[1];
7577841947SLeila Ghaffari   setup_context->dc_axis[2] = dc_axis[2];
7677841947SLeila Ghaffari 
7777841947SLeila Ghaffari   PetscFunctionReturn(0);
7877841947SLeila Ghaffari }
7977841947SLeila Ghaffari 
80d0b732dbSLeila Ghaffari PetscErrorCode SetupContext_DENSITY_CURRENT(Ceed ceed, CeedData ceed_data,
812fe7aee7SLeila Ghaffari     AppCtx app_ctx, SetupContext setup_ctx, Physics phys) {
82d0b732dbSLeila Ghaffari   PetscFunctionBeginUser;
8388b783a1SJames Wright   PetscInt ierr = SetupContext_NEWTONIAN_IG(ceed, ceed_data, app_ctx, setup_ctx,
8488b783a1SJames Wright                   phys);
8588b783a1SJames Wright   CHKERRQ(ierr);
8677841947SLeila Ghaffari   PetscFunctionReturn(0);
8777841947SLeila Ghaffari }
8877841947SLeila Ghaffari 
8977841947SLeila Ghaffari PetscErrorCode PRINT_DENSITY_CURRENT(Physics phys, SetupContext setup_ctx,
9077841947SLeila Ghaffari                                      AppCtx app_ctx) {
9177841947SLeila Ghaffari   MPI_Comm       comm = PETSC_COMM_WORLD;
9277841947SLeila Ghaffari   PetscErrorCode ierr;
9377841947SLeila Ghaffari   PetscFunctionBeginUser;
9477841947SLeila Ghaffari 
9577841947SLeila Ghaffari   ierr = PetscPrintf(comm,
9677841947SLeila Ghaffari                      "  Problem:\n"
9777841947SLeila Ghaffari                      "    Problem Name                       : %s\n"
9877841947SLeila Ghaffari                      "    Stabilization                      : %s\n",
9977841947SLeila Ghaffari                      app_ctx->problem_name, StabilizationTypes[phys->stab]);
10077841947SLeila Ghaffari   CHKERRQ(ierr);
10177841947SLeila Ghaffari 
10277841947SLeila Ghaffari   PetscFunctionReturn(0);
10377841947SLeila Ghaffari }
104