1ccaff030SJeremy L Thompson // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2ccaff030SJeremy L Thompson // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3ccaff030SJeremy L Thompson // reserved. See files LICENSE and NOTICE for details. 4ccaff030SJeremy L Thompson // 5ccaff030SJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software 6ccaff030SJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral 7ccaff030SJeremy L Thompson // element discretizations for exascale applications. For more information and 8ccaff030SJeremy L Thompson // source code availability see http://github.com/ceed. 9ccaff030SJeremy L Thompson // 10ccaff030SJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11ccaff030SJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office 12ccaff030SJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for 13ccaff030SJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including 14ccaff030SJeremy L Thompson // software, applications, hardware, advanced system engineering and early 15ccaff030SJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative. 16ccaff030SJeremy L Thompson 17ccaff030SJeremy L Thompson // libCEED + PETSc Example: Navier-Stokes 18ccaff030SJeremy L Thompson // 19ccaff030SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve a 20ccaff030SJeremy L Thompson // Navier-Stokes problem. 21ccaff030SJeremy L Thompson // 22ccaff030SJeremy L Thompson // The code is intentionally "raw", using only low-level communication 23ccaff030SJeremy L Thompson // primitives. 24ccaff030SJeremy L Thompson // 25ccaff030SJeremy L Thompson // Build with: 26ccaff030SJeremy L Thompson // 27ccaff030SJeremy L Thompson // make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes 28ccaff030SJeremy L Thompson // 29ccaff030SJeremy L Thompson // Sample runs: 30ccaff030SJeremy L Thompson // 3184d34d69SLeila Ghaffari // ./navierstokes -ceed /cpu/self -problem density_current -degree 1 3228688798Sjeremylt // ./navierstokes -ceed /gpu/cuda -problem advection -degree 1 33ccaff030SJeremy L Thompson // 34*dc8efd83SLeila Ghaffari //TESTARGS(name="dc_explicit") -ceed {ceed_resource} -test -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-explicit.bin 35*dc8efd83SLeila Ghaffari //TESTARGS(name="dc_implicit_stab_none") -ceed {ceed_resource} -test -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-implicit-stab-none.bin 36*dc8efd83SLeila Ghaffari //TESTARGS(name="adv_rotation_explicit_strong") -ceed {ceed_resource} -test -problem advection -strong_form 1 -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-explicit-strong.bin 37*dc8efd83SLeila Ghaffari //TESTARGS(name="adv_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -problem advection -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-implicit-stab-supg.bin 38*dc8efd83SLeila Ghaffari //TESTARGS(name="adv_translation_implicit_stab_su") -ceed {ceed_resource} -test -problem advection -CtauS .3 -stab su -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -problem_advection_wind translation -problem_advection_wind_translation .53,-1.33,-2.65 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-translation-implicit-stab-su.bin 39*dc8efd83SLeila Ghaffari //TESTARGS(name="adv2d_rotation_explicit_strong") -ceed {ceed_resource} -test -problem advection2d -strong_form 1 -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-explicit-strong.bin 40*dc8efd83SLeila Ghaffari //TESTARGS(name="adv2d_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -problem advection2d -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-implicit-stab-supg.bin 41*dc8efd83SLeila Ghaffari //TESTARGS(name="adv2d_translation_implicit_stab_su") -ceed {ceed_resource} -test -problem advection2d -CtauS .3 -stab su -degree 3 -dm_plex_box_faces 1,1,2 -units_kilogram 1e-9 -lx 125 -ly 125 -lz 250 -center 62.5,62.5,187.5 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -problem_advection_wind translation -problem_advection_wind_translation .53,-1.33,0 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-translation-implicit-stab-su.bin 42ccaff030SJeremy L Thompson 43ccaff030SJeremy L Thompson /// @file 44ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc 45ccaff030SJeremy L Thompson 46ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n"; 47ccaff030SJeremy L Thompson 48ccaff030SJeremy L Thompson #include <ceed.h> 493d576824SJeremy L Thompson #include <petscdmplex.h> 50ccaff030SJeremy L Thompson #include <petscsys.h> 513d576824SJeremy L Thompson #include <petscts.h> 523d576824SJeremy L Thompson #include <stdbool.h> 53ccaff030SJeremy L Thompson #include "advection.h" 54ccaff030SJeremy L Thompson #include "advection2d.h" 553d576824SJeremy L Thompson #include "common.h" 56ccaff030SJeremy L Thompson #include "densitycurrent.h" 573d576824SJeremy L Thompson #include "setup-boundary.h" 58ccaff030SJeremy L Thompson 5984d34d69SLeila Ghaffari #if PETSC_VERSION_LT(3,14,0) 6084d34d69SLeila Ghaffari # define DMPlexGetClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexGetClosureIndices(a,b,c,d,f,g,i) 6184d34d69SLeila Ghaffari # define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i) 6284d34d69SLeila Ghaffari #endif 6384d34d69SLeila Ghaffari 643ab4fca6SValeria Barra #if PETSC_VERSION_LT(3,14,0) 653ab4fca6SValeria Barra # define DMAddBoundary(a,b,c,d,e,f,g,h,i,j,k,l) DMAddBoundary(a,b,c,d,e,f,g,h,j,k,l) 663ab4fca6SValeria Barra #endif 673ab4fca6SValeria Barra 6884d34d69SLeila Ghaffari // MemType Options 6984d34d69SLeila Ghaffari static const char *const memTypes[] = { 7084d34d69SLeila Ghaffari "host", 7184d34d69SLeila Ghaffari "device", 7284d34d69SLeila Ghaffari "memType", "CEED_MEM_", NULL 7384d34d69SLeila Ghaffari }; 7484d34d69SLeila Ghaffari 75ccaff030SJeremy L Thompson // Problem Options 76ccaff030SJeremy L Thompson typedef enum { 77ccaff030SJeremy L Thompson NS_DENSITY_CURRENT = 0, 78ccaff030SJeremy L Thompson NS_ADVECTION = 1, 79ccaff030SJeremy L Thompson NS_ADVECTION2D = 2, 80ccaff030SJeremy L Thompson } problemType; 81ccaff030SJeremy L Thompson static const char *const problemTypes[] = { 82ccaff030SJeremy L Thompson "density_current", 83ccaff030SJeremy L Thompson "advection", 84ccaff030SJeremy L Thompson "advection2d", 8584d34d69SLeila Ghaffari "problemType", "NS_", NULL 86ccaff030SJeremy L Thompson }; 87ccaff030SJeremy L Thompson 881184866aSLeila Ghaffari // Wind Options for Advection 891184866aSLeila Ghaffari typedef enum { 901184866aSLeila Ghaffari ADVECTION_WIND_ROTATION = 0, 911184866aSLeila Ghaffari ADVECTION_WIND_TRANSLATION = 1, 921184866aSLeila Ghaffari } WindType; 931184866aSLeila Ghaffari static const char *const WindTypes[] = { 941184866aSLeila Ghaffari "rotation", 951184866aSLeila Ghaffari "translation", 961184866aSLeila Ghaffari "WindType", "ADVECTION_WIND_", NULL 971184866aSLeila Ghaffari }; 981184866aSLeila Ghaffari 99ccaff030SJeremy L Thompson typedef enum { 100ccaff030SJeremy L Thompson STAB_NONE = 0, 101ccaff030SJeremy L Thompson STAB_SU = 1, // Streamline Upwind 102ccaff030SJeremy L Thompson STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin 103ccaff030SJeremy L Thompson } StabilizationType; 104ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = { 10584d34d69SLeila Ghaffari "none", 106ccaff030SJeremy L Thompson "SU", 107ccaff030SJeremy L Thompson "SUPG", 108ccaff030SJeremy L Thompson "StabilizationType", "STAB_", NULL 109ccaff030SJeremy L Thompson }; 110ccaff030SJeremy L Thompson 111ccaff030SJeremy L Thompson // Problem specific data 112ccaff030SJeremy L Thompson typedef struct { 1138b982baeSLeila Ghaffari CeedInt dim, qdatasizeVol, qdatasizeSur; 114ebb4b9bdSLeila Ghaffari CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applyVol_ifunction, 115ebb4b9bdSLeila Ghaffari applySur; 116ccaff030SJeremy L Thompson PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, 117ccaff030SJeremy L Thompson PetscScalar[], void *); 1187659d40cSLeila Ghaffari const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc, 1197659d40cSLeila Ghaffari *applyVol_ifunction_loc, *applySur_loc; 120ccaff030SJeremy L Thompson const bool non_zero_time; 121ccaff030SJeremy L Thompson } problemData; 122ccaff030SJeremy L Thompson 123ccaff030SJeremy L Thompson problemData problemOptions[] = { 124ccaff030SJeremy L Thompson [NS_DENSITY_CURRENT] = { 125ccaff030SJeremy L Thompson .dim = 3, 126ea6e0f84SLeila Ghaffari .qdatasizeVol = 10, 1278b982baeSLeila Ghaffari .qdatasizeSur = 4, 128b0137797SLeila Ghaffari .setupVol = Setup, 129b0137797SLeila Ghaffari .setupVol_loc = Setup_loc, 130356fbf4bSLeila Ghaffari .setupSur = SetupBoundary, 131356fbf4bSLeila Ghaffari .setupSur_loc = SetupBoundary_loc, 132ccaff030SJeremy L Thompson .ics = ICsDC, 133ccaff030SJeremy L Thompson .ics_loc = ICsDC_loc, 134c96c872fSLeila Ghaffari .applyVol_rhs = DC, 135c96c872fSLeila Ghaffari .applyVol_rhs_loc = DC_loc, 136c96c872fSLeila Ghaffari .applyVol_ifunction = IFunction_DC, 137c96c872fSLeila Ghaffari .applyVol_ifunction_loc = IFunction_DC_loc, 138ccaff030SJeremy L Thompson .bc = Exact_DC, 13984d34d69SLeila Ghaffari .non_zero_time = PETSC_FALSE, 140ccaff030SJeremy L Thompson }, 141ccaff030SJeremy L Thompson [NS_ADVECTION] = { 142ccaff030SJeremy L Thompson .dim = 3, 143ea6e0f84SLeila Ghaffari .qdatasizeVol = 10, 1448b982baeSLeila Ghaffari .qdatasizeSur = 4, 145b0137797SLeila Ghaffari .setupVol = Setup, 146b0137797SLeila Ghaffari .setupVol_loc = Setup_loc, 147356fbf4bSLeila Ghaffari .setupSur = SetupBoundary, 148356fbf4bSLeila Ghaffari .setupSur_loc = SetupBoundary_loc, 149ccaff030SJeremy L Thompson .ics = ICsAdvection, 150ccaff030SJeremy L Thompson .ics_loc = ICsAdvection_loc, 151c96c872fSLeila Ghaffari .applyVol_rhs = Advection, 152c96c872fSLeila Ghaffari .applyVol_rhs_loc = Advection_loc, 153c96c872fSLeila Ghaffari .applyVol_ifunction = IFunction_Advection, 154c96c872fSLeila Ghaffari .applyVol_ifunction_loc = IFunction_Advection_loc, 1557659d40cSLeila Ghaffari .applySur = Advection_Sur, 1567659d40cSLeila Ghaffari .applySur_loc = Advection_Sur_loc, 157ccaff030SJeremy L Thompson .bc = Exact_Advection, 15884d34d69SLeila Ghaffari .non_zero_time = PETSC_FALSE, 159ccaff030SJeremy L Thompson }, 160ccaff030SJeremy L Thompson [NS_ADVECTION2D] = { 161ccaff030SJeremy L Thompson .dim = 2, 162ea6e0f84SLeila Ghaffari .qdatasizeVol = 5, 1638b982baeSLeila Ghaffari .qdatasizeSur = 3, 164c96c872fSLeila Ghaffari .setupVol = Setup2d, 165c96c872fSLeila Ghaffari .setupVol_loc = Setup2d_loc, 166b0137797SLeila Ghaffari .setupSur = SetupBoundary2d, 167b0137797SLeila Ghaffari .setupSur_loc = SetupBoundary2d_loc, 168ccaff030SJeremy L Thompson .ics = ICsAdvection2d, 169ccaff030SJeremy L Thompson .ics_loc = ICsAdvection2d_loc, 170c96c872fSLeila Ghaffari .applyVol_rhs = Advection2d, 171c96c872fSLeila Ghaffari .applyVol_rhs_loc = Advection2d_loc, 172c96c872fSLeila Ghaffari .applyVol_ifunction = IFunction_Advection2d, 173c96c872fSLeila Ghaffari .applyVol_ifunction_loc = IFunction_Advection2d_loc, 1747659d40cSLeila Ghaffari .applySur = Advection2d_Sur, 1757659d40cSLeila Ghaffari .applySur_loc = Advection2d_Sur_loc, 176ccaff030SJeremy L Thompson .bc = Exact_Advection2d, 17784d34d69SLeila Ghaffari .non_zero_time = PETSC_TRUE, 178ccaff030SJeremy L Thompson }, 179ccaff030SJeremy L Thompson }; 180ccaff030SJeremy L Thompson 181ccaff030SJeremy L Thompson // PETSc user data 182ccaff030SJeremy L Thompson typedef struct User_ *User; 183ccaff030SJeremy L Thompson typedef struct Units_ *Units; 184ccaff030SJeremy L Thompson 185ccaff030SJeremy L Thompson struct User_ { 186ccaff030SJeremy L Thompson MPI_Comm comm; 187ccaff030SJeremy L Thompson PetscInt outputfreq; 188ccaff030SJeremy L Thompson DM dm; 189ccaff030SJeremy L Thompson DM dmviz; 190ccaff030SJeremy L Thompson Mat interpviz; 191ccaff030SJeremy L Thompson Ceed ceed; 192ccaff030SJeremy L Thompson Units units; 193ccaff030SJeremy L Thompson CeedVector qceed, qdotceed, gceed; 1941e150236SLeila Ghaffari CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction; 195ccaff030SJeremy L Thompson Vec M; 196d99129b9SLeila Ghaffari char outputdir[PETSC_MAX_PATH_LEN]; 197ccaff030SJeremy L Thompson PetscInt contsteps; 198ccaff030SJeremy L Thompson }; 199ccaff030SJeremy L Thompson 200ccaff030SJeremy L Thompson struct Units_ { 201ccaff030SJeremy L Thompson // fundamental units 202ccaff030SJeremy L Thompson PetscScalar meter; 203ccaff030SJeremy L Thompson PetscScalar kilogram; 204ccaff030SJeremy L Thompson PetscScalar second; 205ccaff030SJeremy L Thompson PetscScalar Kelvin; 206ccaff030SJeremy L Thompson // derived units 207ccaff030SJeremy L Thompson PetscScalar Pascal; 208ccaff030SJeremy L Thompson PetscScalar JperkgK; 209ccaff030SJeremy L Thompson PetscScalar mpersquareds; 210ccaff030SJeremy L Thompson PetscScalar WpermK; 211ccaff030SJeremy L Thompson PetscScalar kgpercubicm; 212ccaff030SJeremy L Thompson PetscScalar kgpersquaredms; 213ccaff030SJeremy L Thompson PetscScalar Joulepercubicm; 21416c0476cSLeila Ghaffari PetscScalar Joule; 215ccaff030SJeremy L Thompson }; 216ccaff030SJeremy L Thompson 217ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC; 218ccaff030SJeremy L Thompson struct SimpleBC_ { 2197659d40cSLeila Ghaffari PetscInt nwall, nslip[3]; 2207659d40cSLeila Ghaffari PetscInt walls[6], slips[3][6]; 22184d34d69SLeila Ghaffari PetscBool userbc; 222ccaff030SJeremy L Thompson }; 223ccaff030SJeremy L Thompson 224ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1). 225ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) { 226ccaff030SJeremy L Thompson return i >= 0 ? i : -(i+1); 227ccaff030SJeremy L Thompson } 228ccaff030SJeremy L Thompson 229ccaff030SJeremy L Thompson // Utility function to create local CEED restriction 230ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P, 23184d34d69SLeila Ghaffari CeedInt height, DMLabel domainLabel, CeedInt value, 232ccaff030SJeremy L Thompson CeedElemRestriction *Erestrict) { 233ccaff030SJeremy L Thompson 234ccaff030SJeremy L Thompson PetscSection section; 2351184866aSLeila Ghaffari PetscInt p, Nelem, Ndof, *erestrict, eoffset, nfields, dim, depth; 2360c6c0b13SLeila Ghaffari DMLabel depthLabel; 2370c6c0b13SLeila Ghaffari IS depthIS, iterIS; 23884d34d69SLeila Ghaffari Vec Uloc; 2390c6c0b13SLeila Ghaffari const PetscInt *iterIndices; 240ccaff030SJeremy L Thompson PetscErrorCode ierr; 241ccaff030SJeremy L Thompson 242ccaff030SJeremy L Thompson PetscFunctionBeginUser; 243ccaff030SJeremy L Thompson ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 244da51bdd9SLeila Ghaffari dim -= height; 245ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 246ccaff030SJeremy L Thompson ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr); 247ccaff030SJeremy L Thompson PetscInt ncomp[nfields], fieldoff[nfields+1]; 248ccaff030SJeremy L Thompson fieldoff[0] = 0; 249ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 250ccaff030SJeremy L Thompson ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr); 251ccaff030SJeremy L Thompson fieldoff[f+1] = fieldoff[f] + ncomp[f]; 252ccaff030SJeremy L Thompson } 253ccaff030SJeremy L Thompson 2540c6c0b13SLeila Ghaffari ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr); 2550c6c0b13SLeila Ghaffari ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr); 2560c6c0b13SLeila Ghaffari ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr); 2570c6c0b13SLeila Ghaffari if (domainLabel) { 2580c6c0b13SLeila Ghaffari IS domainIS; 2590c6c0b13SLeila Ghaffari ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr); 2601119eeeeSJed Brown if (domainIS) { // domainIS is non-empty 2610c6c0b13SLeila Ghaffari ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr); 2620c6c0b13SLeila Ghaffari ierr = ISDestroy(&domainIS); CHKERRQ(ierr); 2631119eeeeSJed Brown } else { // domainIS is NULL (empty) 2641119eeeeSJed Brown iterIS = NULL; 2651119eeeeSJed Brown } 2660c6c0b13SLeila Ghaffari ierr = ISDestroy(&depthIS); CHKERRQ(ierr); 2670c6c0b13SLeila Ghaffari } else { 2680c6c0b13SLeila Ghaffari iterIS = depthIS; 2690c6c0b13SLeila Ghaffari } 2701119eeeeSJed Brown if (iterIS) { 2710c6c0b13SLeila Ghaffari ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr); 2720c6c0b13SLeila Ghaffari ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr); 2731119eeeeSJed Brown } else { 2741119eeeeSJed Brown Nelem = 0; 2751119eeeeSJed Brown iterIndices = NULL; 2761119eeeeSJed Brown } 277ccaff030SJeremy L Thompson ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr); 2780c6c0b13SLeila Ghaffari for (p=0,eoffset=0; p<Nelem; p++) { 2790c6c0b13SLeila Ghaffari PetscInt c = iterIndices[p]; 280ccaff030SJeremy L Thompson PetscInt numindices, *indices, nnodes; 28184d34d69SLeila Ghaffari ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE, 28284d34d69SLeila Ghaffari &numindices, &indices, NULL, NULL); 28384d34d69SLeila Ghaffari CHKERRQ(ierr); 28432b5ec5fSJed Brown bool flip = false; 28532b5ec5fSJed Brown if (height > 0) { 28632b5ec5fSJed Brown PetscInt numCells, numFaces, start = -1; 28732b5ec5fSJed Brown const PetscInt *orients, *faces, *cells; 28832b5ec5fSJed Brown ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr); 28932b5ec5fSJed Brown ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr); 290ebb4b9bdSLeila Ghaffari if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 291f259b054Svaleriabarra "Expected one cell in support of exterior face, but got %D cells", 292f259b054Svaleriabarra numCells); 29332b5ec5fSJed Brown ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr); 29432b5ec5fSJed Brown ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr); 29532b5ec5fSJed Brown for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;} 296ebb4b9bdSLeila Ghaffari if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, 297f259b054Svaleriabarra "Could not find face %D in cone of its support", 298f259b054Svaleriabarra c); 29932b5ec5fSJed Brown ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr); 30032b5ec5fSJed Brown if (orients[start] < 0) flip = true; 30132b5ec5fSJed Brown } 30284d34d69SLeila Ghaffari if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF, 30384d34d69SLeila Ghaffari PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D", 30484d34d69SLeila Ghaffari c); 305ccaff030SJeremy L Thompson nnodes = numindices / fieldoff[nfields]; 306ccaff030SJeremy L Thompson for (PetscInt i=0; i<nnodes; i++) { 30732b5ec5fSJed Brown PetscInt ii = i; 30832b5ec5fSJed Brown if (flip) { 30932b5ec5fSJed Brown if (P == nnodes) ii = nnodes - 1 - i; 31032b5ec5fSJed Brown else if (P*P == nnodes) { 31132b5ec5fSJed Brown PetscInt row = i / P, col = i % P; 31232b5ec5fSJed Brown ii = row + col * P; 313ebb4b9bdSLeila Ghaffari } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, 314f259b054Svaleriabarra "No support for flipping point with %D nodes != P (%D) or P^2", 315f259b054Svaleriabarra nnodes, P); 31632b5ec5fSJed Brown } 317ccaff030SJeremy L Thompson // Check that indices are blocked by node and thus can be coalesced as a single field with 318ccaff030SJeremy L Thompson // fieldoff[nfields] = sum(ncomp) components. 319ccaff030SJeremy L Thompson for (PetscInt f=0; f<nfields; f++) { 320ccaff030SJeremy L Thompson for (PetscInt j=0; j<ncomp[f]; j++) { 32132b5ec5fSJed Brown if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j]) 32232b5ec5fSJed Brown != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j) 323ccaff030SJeremy L Thompson SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 324ccaff030SJeremy L Thompson "Cell %D closure indices not interlaced for node %D field %D component %D", 32532b5ec5fSJed Brown c, ii, f, j); 326ccaff030SJeremy L Thompson } 327ccaff030SJeremy L Thompson } 328ccaff030SJeremy L Thompson // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode. 32932b5ec5fSJed Brown PetscInt loc = Involute(indices[ii*ncomp[0]]); 3306f55dfd5Svaleriabarra erestrict[eoffset++] = loc; 331ccaff030SJeremy L Thompson } 33284d34d69SLeila Ghaffari ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE, 33384d34d69SLeila Ghaffari &numindices, &indices, NULL, NULL); 33484d34d69SLeila Ghaffari CHKERRQ(ierr); 335ccaff030SJeremy L Thompson } 3360c6c0b13SLeila Ghaffari if (eoffset != Nelem*PetscPowInt(P, dim)) 3370c6c0b13SLeila Ghaffari SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB, 3380c6c0b13SLeila Ghaffari "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem, 339ccaff030SJeremy L Thompson PetscPowInt(P, dim),eoffset); 3401119eeeeSJed Brown if (iterIS) { 3410c6c0b13SLeila Ghaffari ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr); 3421119eeeeSJed Brown } 3430c6c0b13SLeila Ghaffari ierr = ISDestroy(&iterIS); CHKERRQ(ierr); 3440c6c0b13SLeila Ghaffari 345ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr); 346ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr); 347ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr); 3486f55dfd5Svaleriabarra CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields], 3496f55dfd5Svaleriabarra 1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict, 3506f55dfd5Svaleriabarra Erestrict); 351ccaff030SJeremy L Thompson ierr = PetscFree(erestrict); CHKERRQ(ierr); 352ccaff030SJeremy L Thompson PetscFunctionReturn(0); 353ccaff030SJeremy L Thompson } 354ccaff030SJeremy L Thompson 355c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain 3561e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, 3571e150236SLeila Ghaffari DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize, 3581e150236SLeila Ghaffari CeedElemRestriction *restrictq, CeedElemRestriction *restrictx, 3591e150236SLeila Ghaffari CeedElemRestriction *restrictqdi) { 360c96c872fSLeila Ghaffari 361c96c872fSLeila Ghaffari DM dmcoord; 3621e150236SLeila Ghaffari CeedInt dim, localNelem; 3631e150236SLeila Ghaffari CeedInt Qdim; 364c96c872fSLeila Ghaffari PetscErrorCode ierr; 365c96c872fSLeila Ghaffari 366c96c872fSLeila Ghaffari PetscFunctionBeginUser; 3671e150236SLeila Ghaffari ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 3681e150236SLeila Ghaffari dim -= height; 3691e150236SLeila Ghaffari Qdim = CeedIntPow(Q, dim); 370c96c872fSLeila Ghaffari ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 371c96c872fSLeila Ghaffari ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 372c96c872fSLeila Ghaffari CHKERRQ(ierr); 373ebb4b9bdSLeila Ghaffari ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, 374ebb4b9bdSLeila Ghaffari restrictq); 375c96c872fSLeila Ghaffari CHKERRQ(ierr); 376ebb4b9bdSLeila Ghaffari ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, 377ebb4b9bdSLeila Ghaffari restrictx); 378c96c872fSLeila Ghaffari CHKERRQ(ierr); 379c96c872fSLeila Ghaffari CeedElemRestrictionGetNumElements(*restrictq, &localNelem); 380c96c872fSLeila Ghaffari CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim, 381c96c872fSLeila Ghaffari qdatasize, qdatasize*localNelem*Qdim, 382c96c872fSLeila Ghaffari CEED_STRIDES_BACKEND, restrictqdi); 383c96c872fSLeila Ghaffari PetscFunctionReturn(0); 384c96c872fSLeila Ghaffari } 385c96c872fSLeila Ghaffari 3861e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain 3877659d40cSLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, 3887659d40cSLeila Ghaffari WindType wind_type, CeedOperator op_applyVol, CeedQFunction qf_applySur, 3897659d40cSLeila Ghaffari CeedQFunction qf_setupSur, CeedInt height, CeedInt numP_Sur, CeedInt numQ_Sur, 3907659d40cSLeila Ghaffari CeedInt qdatasizeSur, CeedInt NqptsSur, CeedBasis basisxSur, 3917659d40cSLeila Ghaffari CeedBasis basisqSur, CeedOperator *op_apply) { 392ca3ac6ddSLeila Ghaffari 3937659d40cSLeila Ghaffari CeedInt dim, nFace; 3947659d40cSLeila Ghaffari PetscInt lsize, localNelemSur[6]; 3951e150236SLeila Ghaffari Vec Xloc; 3967659d40cSLeila Ghaffari CeedVector xcorners, qdataSur[6]; 3977659d40cSLeila Ghaffari CeedOperator op_setupSur[6], op_applySur[6]; 3987659d40cSLeila Ghaffari CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6]; 399ca3ac6ddSLeila Ghaffari DMLabel domainLabel; 4001e150236SLeila Ghaffari PetscScalar *x; 401ca3ac6ddSLeila Ghaffari PetscErrorCode ierr; 402ca3ac6ddSLeila Ghaffari 403ca3ac6ddSLeila Ghaffari PetscFunctionBeginUser; 404ca3ac6ddSLeila Ghaffari // Composite Operaters 405ca3ac6ddSLeila Ghaffari CeedCompositeOperatorCreate(ceed, op_apply); 406ebb4b9bdSLeila Ghaffari CeedCompositeOperatorAddSub(*op_apply, 407ebb4b9bdSLeila Ghaffari op_applyVol); // Apply a Sub-Operator for the volume 408ca3ac6ddSLeila Ghaffari 4097659d40cSLeila Ghaffari if (wind_type == ADVECTION_WIND_TRANSLATION) { 4107659d40cSLeila Ghaffari bc->nwall = 0; 4117659d40cSLeila Ghaffari bc->nslip[0] = bc->nslip[1] = bc->nslip[2] = 0; 4121e150236SLeila Ghaffari ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 4131e150236SLeila Ghaffari ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr); 4141e150236SLeila Ghaffari ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr); 4151e150236SLeila Ghaffari ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr); 4161e150236SLeila Ghaffari CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x); 417ca3ac6ddSLeila Ghaffari ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr); 4187659d40cSLeila Ghaffari ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr); 4197659d40cSLeila Ghaffari if (dim == 2) nFace = 4; 4207659d40cSLeila Ghaffari if (dim == 3) nFace = 6; 4219fe13df9SLeila Ghaffari 4227659d40cSLeila Ghaffari // Create CEED Operator for each boundary face 4237659d40cSLeila Ghaffari for (CeedInt i=0; i<nFace; i++) { 4247659d40cSLeila Ghaffari ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, i+1, numP_Sur, 425f259b054Svaleriabarra numQ_Sur, qdatasizeSur, &restrictqSur[i], 426f259b054Svaleriabarra &restrictxSur[i], &restrictqdiSur[i]); 427f259b054Svaleriabarra CHKERRQ(ierr); 4287659d40cSLeila Ghaffari // Create the CEED vectors that will be needed in Boundary setup 4297659d40cSLeila Ghaffari CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]); 430f259b054Svaleriabarra CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, 431f259b054Svaleriabarra &qdataSur[i]); 4327659d40cSLeila Ghaffari // Create the operator that builds the quadrature data for the Boundary operator 4337659d40cSLeila Ghaffari CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]); 434ebb4b9bdSLeila Ghaffari CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur, 435ebb4b9bdSLeila Ghaffari CEED_VECTOR_ACTIVE); 4367659d40cSLeila Ghaffari CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE, 437ca3ac6ddSLeila Ghaffari basisxSur, CEED_VECTOR_NONE); 4387659d40cSLeila Ghaffari CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i], 439ca3ac6ddSLeila Ghaffari CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 4407659d40cSLeila Ghaffari // Create Boundary operator 4417659d40cSLeila Ghaffari CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]); 442ebb4b9bdSLeila Ghaffari CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur, 443ebb4b9bdSLeila Ghaffari CEED_VECTOR_ACTIVE); 4447659d40cSLeila Ghaffari CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i], 4457659d40cSLeila Ghaffari CEED_BASIS_COLLOCATED, qdataSur[i]); 446f259b054Svaleriabarra CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur, 447f259b054Svaleriabarra xcorners); 448ebb4b9bdSLeila Ghaffari CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur, 449ebb4b9bdSLeila Ghaffari CEED_VECTOR_ACTIVE); 4507659d40cSLeila Ghaffari // Apply CEED operator for Boundary setup 451ebb4b9bdSLeila Ghaffari CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i], 452ebb4b9bdSLeila Ghaffari CEED_REQUEST_IMMEDIATE); 4537659d40cSLeila Ghaffari // Apply Sub-Operator for the Boundary 4547659d40cSLeila Ghaffari CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]); 4559fe13df9SLeila Ghaffari } 4561e150236SLeila Ghaffari CeedVectorDestroy(&xcorners); 457ca3ac6ddSLeila Ghaffari } 458ca3ac6ddSLeila Ghaffari PetscFunctionReturn(0); 459ca3ac6ddSLeila Ghaffari } 460ca3ac6ddSLeila Ghaffari 461ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) { 462ccaff030SJeremy L Thompson PetscErrorCode ierr; 463ccaff030SJeremy L Thompson PetscInt m; 464ccaff030SJeremy L Thompson 465ccaff030SJeremy L Thompson PetscFunctionBeginUser; 466ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr); 467ccaff030SJeremy L Thompson ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr); 468ccaff030SJeremy L Thompson PetscFunctionReturn(0); 469ccaff030SJeremy L Thompson } 470ccaff030SJeremy L Thompson 471ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) { 472ccaff030SJeremy L Thompson PetscErrorCode ierr; 473ccaff030SJeremy L Thompson PetscInt mceed, mpetsc; 474ccaff030SJeremy L Thompson PetscScalar *a; 475ccaff030SJeremy L Thompson 476ccaff030SJeremy L Thompson PetscFunctionBeginUser; 477ccaff030SJeremy L Thompson ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr); 478ccaff030SJeremy L Thompson ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr); 479ccaff030SJeremy L Thompson if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, 48084d34d69SLeila Ghaffari "Cannot place PETSc Vec of length %D in CeedVector of length %D", 48184d34d69SLeila Ghaffari mpetsc, mceed); 482ccaff030SJeremy L Thompson ierr = VecGetArray(p, &a); CHKERRQ(ierr); 483ccaff030SJeremy L Thompson CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a); 484ccaff030SJeremy L Thompson PetscFunctionReturn(0); 485ccaff030SJeremy L Thompson } 486ccaff030SJeremy L Thompson 487ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, 488ccaff030SJeremy L Thompson PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM, 489ccaff030SJeremy L Thompson Vec cellGeomFVM, Vec gradFVM) { 490ccaff030SJeremy L Thompson PetscErrorCode ierr; 491ccaff030SJeremy L Thompson Vec Qbc; 492ccaff030SJeremy L Thompson 493ccaff030SJeremy L Thompson PetscFunctionBegin; 494ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 495ccaff030SJeremy L Thompson ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr); 496ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 497ccaff030SJeremy L Thompson PetscFunctionReturn(0); 498ccaff030SJeremy L Thompson } 499ccaff030SJeremy L Thompson 500ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u) 501ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G 502ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) { 503ccaff030SJeremy L Thompson PetscErrorCode ierr; 504ccaff030SJeremy L Thompson User user = *(User *)userData; 505ccaff030SJeremy L Thompson PetscScalar *q, *g; 506ccaff030SJeremy L Thompson Vec Qloc, Gloc; 507ccaff030SJeremy L Thompson 508ccaff030SJeremy L Thompson // Global-to-local 509ccaff030SJeremy L Thompson PetscFunctionBeginUser; 510ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 511ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 512ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 513ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 514ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 515ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 516ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 517ccaff030SJeremy L Thompson 518ccaff030SJeremy L Thompson // Ceed Vectors 519ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 520ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 521ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q); 522ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 523ccaff030SJeremy L Thompson 524ccaff030SJeremy L Thompson // Apply CEED operator 525ccaff030SJeremy L Thompson CeedOperatorApply(user->op_rhs, user->qceed, user->gceed, 526ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 527ccaff030SJeremy L Thompson 528ccaff030SJeremy L Thompson // Restore vectors 529ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr); 530ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 531ccaff030SJeremy L Thompson 532ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 533ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 534ccaff030SJeremy L Thompson 535ccaff030SJeremy L Thompson // Inverse of the lumped mass matrix 536ccaff030SJeremy L Thompson ierr = VecPointwiseMult(G, G, user->M); // M is Minv 537ccaff030SJeremy L Thompson CHKERRQ(ierr); 538ccaff030SJeremy L Thompson 539ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 540ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 541ccaff030SJeremy L Thompson PetscFunctionReturn(0); 542ccaff030SJeremy L Thompson } 543ccaff030SJeremy L Thompson 544ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G, 545ccaff030SJeremy L Thompson void *userData) { 546ccaff030SJeremy L Thompson PetscErrorCode ierr; 547ccaff030SJeremy L Thompson User user = *(User *)userData; 548ccaff030SJeremy L Thompson const PetscScalar *q, *qdot; 549ccaff030SJeremy L Thompson PetscScalar *g; 550ccaff030SJeremy L Thompson Vec Qloc, Qdotloc, Gloc; 551ccaff030SJeremy L Thompson 552ccaff030SJeremy L Thompson // Global-to-local 553ccaff030SJeremy L Thompson PetscFunctionBeginUser; 554ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 555ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 556ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 557ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 558ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 559ccaff030SJeremy L Thompson ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0, 560ccaff030SJeremy L Thompson NULL, NULL, NULL); CHKERRQ(ierr); 561ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr); 562ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr); 563ccaff030SJeremy L Thompson ierr = VecZeroEntries(Gloc); CHKERRQ(ierr); 564ccaff030SJeremy L Thompson 565ccaff030SJeremy L Thompson // Ceed Vectors 566ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr); 567ccaff030SJeremy L Thompson ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 568ccaff030SJeremy L Thompson ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr); 569ccaff030SJeremy L Thompson CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, 570ccaff030SJeremy L Thompson (PetscScalar *)q); 571ccaff030SJeremy L Thompson CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER, 572ccaff030SJeremy L Thompson (PetscScalar *)qdot); 573ccaff030SJeremy L Thompson CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g); 574ccaff030SJeremy L Thompson 575ccaff030SJeremy L Thompson // Apply CEED operator 576ccaff030SJeremy L Thompson CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed, 577ccaff030SJeremy L Thompson CEED_REQUEST_IMMEDIATE); 578ccaff030SJeremy L Thompson 579ccaff030SJeremy L Thompson // Restore vectors 580ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr); 581ccaff030SJeremy L Thompson ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr); 582ccaff030SJeremy L Thompson ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr); 583ccaff030SJeremy L Thompson 584ccaff030SJeremy L Thompson ierr = VecZeroEntries(G); CHKERRQ(ierr); 585ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr); 586ccaff030SJeremy L Thompson 587ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 588ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr); 589ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr); 590ccaff030SJeremy L Thompson PetscFunctionReturn(0); 591ccaff030SJeremy L Thompson } 592ccaff030SJeremy L Thompson 593ccaff030SJeremy L Thompson // User provided TS Monitor 594ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time, 595ccaff030SJeremy L Thompson Vec Q, void *ctx) { 596ccaff030SJeremy L Thompson User user = ctx; 597ccaff030SJeremy L Thompson Vec Qloc; 598ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 599ccaff030SJeremy L Thompson PetscViewer viewer; 600ccaff030SJeremy L Thompson PetscErrorCode ierr; 601ccaff030SJeremy L Thompson 602ccaff030SJeremy L Thompson // Set up output 603ccaff030SJeremy L Thompson PetscFunctionBeginUser; 604ccaff030SJeremy L Thompson // Print every 'outputfreq' steps 605ccaff030SJeremy L Thompson if (stepno % user->outputfreq != 0) 606ccaff030SJeremy L Thompson PetscFunctionReturn(0); 607ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 608ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr); 609ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 610ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 611ccaff030SJeremy L Thompson 612ccaff030SJeremy L Thompson // Output 613ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu", 614d99129b9SLeila Ghaffari user->outputdir, stepno + user->contsteps); 615ccaff030SJeremy L Thompson CHKERRQ(ierr); 616ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath, 617ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer); CHKERRQ(ierr); 618ccaff030SJeremy L Thompson ierr = VecView(Qloc, viewer); CHKERRQ(ierr); 6199d801c56SJed Brown ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 620ccaff030SJeremy L Thompson if (user->dmviz) { 621ccaff030SJeremy L Thompson Vec Qrefined, Qrefined_loc; 622ccaff030SJeremy L Thompson char filepath_refined[PETSC_MAX_PATH_LEN]; 623ccaff030SJeremy L Thompson PetscViewer viewer_refined; 624ccaff030SJeremy L Thompson 625ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 626ccaff030SJeremy L Thompson ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 627ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined"); 628ccaff030SJeremy L Thompson CHKERRQ(ierr); 629ccaff030SJeremy L Thompson ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr); 630ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr); 631ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc); 632ccaff030SJeremy L Thompson CHKERRQ(ierr); 633ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined, 634ccaff030SJeremy L Thompson "%s/nsrefined-%03D.vtu", 635d99129b9SLeila Ghaffari user->outputdir, stepno + user->contsteps); 636ccaff030SJeremy L Thompson CHKERRQ(ierr); 637ccaff030SJeremy L Thompson ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined), 638ccaff030SJeremy L Thompson filepath_refined, 639ccaff030SJeremy L Thompson FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr); 640ccaff030SJeremy L Thompson ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr); 641ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr); 642ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr); 643ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr); 644ccaff030SJeremy L Thompson } 645ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr); 646ccaff030SJeremy L Thompson 647ccaff030SJeremy L Thompson // Save data in a binary file for continuation of simulations 648ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 649d99129b9SLeila Ghaffari user->outputdir); CHKERRQ(ierr); 650ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 651ccaff030SJeremy L Thompson CHKERRQ(ierr); 652ccaff030SJeremy L Thompson ierr = VecView(Q, viewer); CHKERRQ(ierr); 653ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 654ccaff030SJeremy L Thompson 655ccaff030SJeremy L Thompson // Save time stamp 656ccaff030SJeremy L Thompson // Dimensionalize time back 657ccaff030SJeremy L Thompson time /= user->units->second; 658ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 659d99129b9SLeila Ghaffari user->outputdir); CHKERRQ(ierr); 660ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer); 661ccaff030SJeremy L Thompson CHKERRQ(ierr); 662ccaff030SJeremy L Thompson #if PETSC_VERSION_GE(3,13,0) 663ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL); 664ccaff030SJeremy L Thompson #else 665ccaff030SJeremy L Thompson ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true); 666ccaff030SJeremy L Thompson #endif 667ccaff030SJeremy L Thompson CHKERRQ(ierr); 668ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 669ccaff030SJeremy L Thompson 670ccaff030SJeremy L Thompson PetscFunctionReturn(0); 671ccaff030SJeremy L Thompson } 672ccaff030SJeremy L Thompson 67384d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics, 674ccaff030SJeremy L Thompson CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q, 675777ff853SJeremy L Thompson CeedElemRestriction restrictq, CeedQFunctionContext ctxSetup, CeedScalar time) { 676ccaff030SJeremy L Thompson PetscErrorCode ierr; 677ccaff030SJeremy L Thompson CeedVector multlvec; 678ccaff030SJeremy L Thompson Vec Multiplicity, MultiplicityLoc; 679ccaff030SJeremy L Thompson 680777ff853SJeremy L Thompson SetupContext ctxSetupData; 681777ff853SJeremy L Thompson CeedQFunctionContextGetData(ctxSetup, CEED_MEM_HOST, (void **)&ctxSetupData); 682777ff853SJeremy L Thompson ctxSetupData->time = time; 683777ff853SJeremy L Thompson CeedQFunctionContextRestoreData(ctxSetup, (void **)&ctxSetupData); 684777ff853SJeremy L Thompson 685ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 686ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 687ccaff030SJeremy L Thompson CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE); 688ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 689ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr); 690ccaff030SJeremy L Thompson 691ccaff030SJeremy L Thompson // Fix multiplicity for output of ICs 692ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 693ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL); 694ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr); 695ccaff030SJeremy L Thompson CeedElemRestrictionGetMultiplicity(restrictq, multlvec); 696ccaff030SJeremy L Thompson CeedVectorDestroy(&multlvec); 697ccaff030SJeremy L Thompson ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 698ccaff030SJeremy L Thompson ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr); 699ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity); 700ccaff030SJeremy L Thompson CHKERRQ(ierr); 701ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr); 702ccaff030SJeremy L Thompson ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr); 703ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr); 704ccaff030SJeremy L Thompson ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr); 705ccaff030SJeremy L Thompson 706ccaff030SJeremy L Thompson PetscFunctionReturn(0); 707ccaff030SJeremy L Thompson } 708ccaff030SJeremy L Thompson 709ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, 710ccaff030SJeremy L Thompson CeedElemRestriction restrictq, CeedBasis basisq, 711ccaff030SJeremy L Thompson CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) { 712ccaff030SJeremy L Thompson PetscErrorCode ierr; 713ccaff030SJeremy L Thompson CeedQFunction qf_mass; 714ccaff030SJeremy L Thompson CeedOperator op_mass; 715ccaff030SJeremy L Thompson CeedVector mceed; 716ccaff030SJeremy L Thompson Vec Mloc; 717ccaff030SJeremy L Thompson CeedInt ncompq, qdatasize; 718ccaff030SJeremy L Thompson 719ccaff030SJeremy L Thompson PetscFunctionBeginUser; 720ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictq, &ncompq); 721ccaff030SJeremy L Thompson CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize); 722ccaff030SJeremy L Thompson // Create the Q-function that defines the action of the mass operator 723ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass); 724ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP); 725ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE); 726ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP); 727ccaff030SJeremy L Thompson 728ccaff030SJeremy L Thompson // Create the mass operator 729ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 730ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 731ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "qdata", restrictqdi, 732ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, qdata); 733ccaff030SJeremy L Thompson CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 734ccaff030SJeremy L Thompson 735ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr); 736ccaff030SJeremy L Thompson ierr = VecZeroEntries(Mloc); CHKERRQ(ierr); 737ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &mceed, NULL); 738ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr); 739ccaff030SJeremy L Thompson 740ccaff030SJeremy L Thompson { 741ccaff030SJeremy L Thompson // Compute a lumped mass matrix 742ccaff030SJeremy L Thompson CeedVector onesvec; 743ccaff030SJeremy L Thompson CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL); 744ccaff030SJeremy L Thompson CeedVectorSetValue(onesvec, 1.0); 745ccaff030SJeremy L Thompson CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE); 746ccaff030SJeremy L Thompson CeedVectorDestroy(&onesvec); 747ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_mass); 748ccaff030SJeremy L Thompson CeedVectorDestroy(&mceed); 749ccaff030SJeremy L Thompson } 750ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_mass); 751ccaff030SJeremy L Thompson 752ccaff030SJeremy L Thompson ierr = VecZeroEntries(M); CHKERRQ(ierr); 753ccaff030SJeremy L Thompson ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr); 754ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr); 755ccaff030SJeremy L Thompson 756ccaff030SJeremy L Thompson // Invert diagonally lumped mass vector for RHS function 757ccaff030SJeremy L Thompson ierr = VecReciprocal(M); CHKERRQ(ierr); 758ccaff030SJeremy L Thompson PetscFunctionReturn(0); 759ccaff030SJeremy L Thompson } 760ccaff030SJeremy L Thompson 76184d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree, 762777ff853SJeremy L Thompson SimpleBC bc, void *ctxSetupData) { 763ccaff030SJeremy L Thompson PetscErrorCode ierr; 764ccaff030SJeremy L Thompson 765ccaff030SJeremy L Thompson PetscFunctionBeginUser; 766ccaff030SJeremy L Thompson { 767ccaff030SJeremy L Thompson // Configure the finite element space and boundary conditions 768ccaff030SJeremy L Thompson PetscFE fe; 769ccaff030SJeremy L Thompson PetscInt ncompq = 5; 770ff6701fcSJed Brown ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq, 771ff6701fcSJed Brown PETSC_FALSE, degree, PETSC_DECIDE, 77232ed2d11SJed Brown &fe); CHKERRQ(ierr); 773ccaff030SJeremy L Thompson ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr); 774ccaff030SJeremy L Thompson ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr); 775ccaff030SJeremy L Thompson ierr = DMCreateDS(dm); CHKERRQ(ierr); 77607af6069Svaleriabarra { 77707af6069Svaleriabarra PetscInt comps[1] = {1}; 77807af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0, 7793ab4fca6SValeria Barra 1, comps, (void(*)(void))NULL, NULL, bc->nslip[0], 780777ff853SJeremy L Thompson bc->slips[0], ctxSetupData); CHKERRQ(ierr); 78107af6069Svaleriabarra comps[0] = 2; 78207af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0, 7833ab4fca6SValeria Barra 1, comps, (void(*)(void))NULL, NULL, bc->nslip[1], 784777ff853SJeremy L Thompson bc->slips[1], ctxSetupData); CHKERRQ(ierr); 78507af6069Svaleriabarra comps[0] = 3; 78607af6069Svaleriabarra ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0, 7873ab4fca6SValeria Barra 1, comps, (void(*)(void))NULL, NULL, bc->nslip[2], 788777ff853SJeremy L Thompson bc->slips[2], ctxSetupData); CHKERRQ(ierr); 78907af6069Svaleriabarra } 79084d34d69SLeila Ghaffari if (bc->userbc == PETSC_TRUE) { 79184d34d69SLeila Ghaffari for (PetscInt c = 0; c < 3; c++) { 79284d34d69SLeila Ghaffari for (PetscInt s = 0; s < bc->nslip[c]; s++) { 79384d34d69SLeila Ghaffari for (PetscInt w = 0; w < bc->nwall; w++) { 79484d34d69SLeila Ghaffari if (bc->slips[c][s] == bc->walls[w]) 79584d34d69SLeila Ghaffari SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, 796f259b054Svaleriabarra "Boundary condition already set on face %D!\n", 797f259b054Svaleriabarra bc->walls[w]); 79884d34d69SLeila Ghaffari 79984d34d69SLeila Ghaffari } 80084d34d69SLeila Ghaffari } 80184d34d69SLeila Ghaffari } 80284d34d69SLeila Ghaffari } 80384d34d69SLeila Ghaffari // Wall boundary conditions are zero energy density and zero flux for 80484d34d69SLeila Ghaffari // velocity in advection/advection2d, and zero velocity and zero flux 80584d34d69SLeila Ghaffari // for mass density and energy density in density_current 80684d34d69SLeila Ghaffari { 80784d34d69SLeila Ghaffari if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) { 80884d34d69SLeila Ghaffari PetscInt comps[1] = {4}; 80984d34d69SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 8103ab4fca6SValeria Barra 1, comps, (void(*)(void))problem->bc, NULL, 811777ff853SJeremy L Thompson bc->nwall, bc->walls, ctxSetupData); CHKERRQ(ierr); 81284d34d69SLeila Ghaffari } else if (problem->bc == Exact_DC) { 81384d34d69SLeila Ghaffari PetscInt comps[3] = {1, 2, 3}; 81484d34d69SLeila Ghaffari ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0, 8153ab4fca6SValeria Barra 3, comps, (void(*)(void))problem->bc, NULL, 816777ff853SJeremy L Thompson bc->nwall, bc->walls, ctxSetupData); CHKERRQ(ierr); 81784d34d69SLeila Ghaffari } else 81884d34d69SLeila Ghaffari SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, 81984d34d69SLeila Ghaffari "Undefined boundary conditions for this problem"); 82084d34d69SLeila Ghaffari } 821ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL); 822ccaff030SJeremy L Thompson CHKERRQ(ierr); 823ccaff030SJeremy L Thompson ierr = PetscFEDestroy(&fe); CHKERRQ(ierr); 824ccaff030SJeremy L Thompson } 825ccaff030SJeremy L Thompson { 826ccaff030SJeremy L Thompson // Empty name for conserved field (because there is only one field) 827ccaff030SJeremy L Thompson PetscSection section; 828ccaff030SJeremy L Thompson ierr = DMGetLocalSection(dm, §ion); CHKERRQ(ierr); 829ccaff030SJeremy L Thompson ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr); 830ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 0, "Density"); 831ccaff030SJeremy L Thompson CHKERRQ(ierr); 832ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX"); 833ccaff030SJeremy L Thompson CHKERRQ(ierr); 834ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY"); 835ccaff030SJeremy L Thompson CHKERRQ(ierr); 836ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ"); 837ccaff030SJeremy L Thompson CHKERRQ(ierr); 838ccaff030SJeremy L Thompson ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity"); 839ccaff030SJeremy L Thompson CHKERRQ(ierr); 840ccaff030SJeremy L Thompson } 841ccaff030SJeremy L Thompson PetscFunctionReturn(0); 842ccaff030SJeremy L Thompson } 843ccaff030SJeremy L Thompson 844ccaff030SJeremy L Thompson int main(int argc, char **argv) { 845ccaff030SJeremy L Thompson PetscInt ierr; 846ccaff030SJeremy L Thompson MPI_Comm comm; 84784d34d69SLeila Ghaffari DM dm, dmcoord, dmviz; 848ccaff030SJeremy L Thompson Mat interpviz; 849ccaff030SJeremy L Thompson TS ts; 850ccaff030SJeremy L Thompson TSAdapt adapt; 851ccaff030SJeremy L Thompson User user; 852ccaff030SJeremy L Thompson Units units; 853ccaff030SJeremy L Thompson char ceedresource[4096] = "/cpu/self"; 85484d34d69SLeila Ghaffari PetscInt localNelemVol, lnodes, gnodes, steps; 855ccaff030SJeremy L Thompson const PetscInt ncompq = 5; 856ccaff030SJeremy L Thompson PetscMPIInt rank; 857ccaff030SJeremy L Thompson PetscScalar ftime; 858ccaff030SJeremy L Thompson Vec Q, Qloc, Xloc; 859ccaff030SJeremy L Thompson Ceed ceed; 86084d34d69SLeila Ghaffari CeedInt numP, numQ; 861cfa64770SLeila Ghaffari CeedVector xcorners, qdata, q0ceed; 86284d34d69SLeila Ghaffari CeedBasis basisx, basisxc, basisq; 86384d34d69SLeila Ghaffari CeedElemRestriction restrictx, restrictq, restrictqdi; 864cfa64770SLeila Ghaffari CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol; 865777ff853SJeremy L Thompson CeedQFunctionContext ctxSetup, ctxNS, ctxAdvection2d, ctxSurface; 866cfa64770SLeila Ghaffari CeedOperator op_setupVol, op_ics; 867ccaff030SJeremy L Thompson CeedScalar Rd; 86884d34d69SLeila Ghaffari CeedMemType memtyperequested; 869ccaff030SJeremy L Thompson PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm, 87016c0476cSLeila Ghaffari kgpersquaredms, Joulepercubicm, Joule; 871ccaff030SJeremy L Thompson problemType problemChoice; 872ccaff030SJeremy L Thompson problemData *problem = NULL; 8731184866aSLeila Ghaffari WindType wind_type; 874ccaff030SJeremy L Thompson StabilizationType stab; 87584d34d69SLeila Ghaffari PetscBool implicit; 876cb3e2689Svaleriabarra PetscInt viz_refine = 0; 877ccaff030SJeremy L Thompson struct SimpleBC_ bc = { 87884d34d69SLeila Ghaffari .nslip = {2, 2, 2}, 87984d34d69SLeila Ghaffari .slips = {{5, 6}, {3, 4}, {1, 2}} 880ccaff030SJeremy L Thompson }; 881ccaff030SJeremy L Thompson double start, cpu_time_used; 882*dc8efd83SLeila Ghaffari // Test variables 883*dc8efd83SLeila Ghaffari PetscBool test; 884*dc8efd83SLeila Ghaffari PetscScalar testtol = 0.; 885*dc8efd83SLeila Ghaffari char filepath[PETSC_MAX_PATH_LEN]; 88684d34d69SLeila Ghaffari // Check PETSc CUDA support 88784d34d69SLeila Ghaffari PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE; 88884d34d69SLeila Ghaffari // *INDENT-OFF* 88984d34d69SLeila Ghaffari #ifdef PETSC_HAVE_CUDA 89084d34d69SLeila Ghaffari petschavecuda = PETSC_TRUE; 89184d34d69SLeila Ghaffari #else 89284d34d69SLeila Ghaffari petschavecuda = PETSC_FALSE; 89384d34d69SLeila Ghaffari #endif 89484d34d69SLeila Ghaffari // *INDENT-ON* 895ccaff030SJeremy L Thompson 896ccaff030SJeremy L Thompson // Create the libCEED contexts 897ccaff030SJeremy L Thompson PetscScalar meter = 1e-2; // 1 meter in scaled length units 898ccaff030SJeremy L Thompson PetscScalar second = 1e-2; // 1 second in scaled time units 899ccaff030SJeremy L Thompson PetscScalar kilogram = 1e-6; // 1 kilogram in scaled mass units 900ccaff030SJeremy L Thompson PetscScalar Kelvin = 1; // 1 Kelvin in scaled temperature units 901ccaff030SJeremy L Thompson CeedScalar theta0 = 300.; // K 902ccaff030SJeremy L Thompson CeedScalar thetaC = -15.; // K 903ccaff030SJeremy L Thompson CeedScalar P0 = 1.e5; // Pa 90416c0476cSLeila Ghaffari CeedScalar E_wind = 1.e6; // J 905ccaff030SJeremy L Thompson CeedScalar N = 0.01; // 1/s 906ccaff030SJeremy L Thompson CeedScalar cv = 717.; // J/(kg K) 907ccaff030SJeremy L Thompson CeedScalar cp = 1004.; // J/(kg K) 908ccaff030SJeremy L Thompson CeedScalar g = 9.81; // m/s^2 909ccaff030SJeremy L Thompson CeedScalar lambda = -2./3.; // - 910ccaff030SJeremy L Thompson CeedScalar mu = 75.; // Pa s, dynamic viscosity 911ccaff030SJeremy L Thompson // mu = 75 is not physical for air, but is good for numerical stability 912ccaff030SJeremy L Thompson CeedScalar k = 0.02638; // W/(m K) 913ccaff030SJeremy L Thompson CeedScalar CtauS = 0.; // dimensionless 914ccaff030SJeremy L Thompson CeedScalar strong_form = 0.; // [0,1] 915ccaff030SJeremy L Thompson PetscScalar lx = 8000.; // m 916ccaff030SJeremy L Thompson PetscScalar ly = 8000.; // m 917ccaff030SJeremy L Thompson PetscScalar lz = 4000.; // m 918ccaff030SJeremy L Thompson CeedScalar rc = 1000.; // m (Radius of bubble) 919ccaff030SJeremy L Thompson PetscScalar resx = 1000.; // m (resolution in x) 920ccaff030SJeremy L Thompson PetscScalar resy = 1000.; // m (resolution in y) 921ccaff030SJeremy L Thompson PetscScalar resz = 1000.; // m (resolution in z) 922ccaff030SJeremy L Thompson PetscInt outputfreq = 10; // - 923ccaff030SJeremy L Thompson PetscInt contsteps = 0; // - 92484d34d69SLeila Ghaffari PetscInt degree = 1; // - 92584d34d69SLeila Ghaffari PetscInt qextra = 2; // - 926ea6e0f84SLeila Ghaffari PetscInt qextraSur = 2; // - 9271184866aSLeila Ghaffari PetscReal center[3], dc_axis[3] = {0, 0, 0}, wind[3] = {1., 0, 0}; 928ccaff030SJeremy L Thompson 929ccaff030SJeremy L Thompson ierr = PetscInitialize(&argc, &argv, NULL, help); 930ccaff030SJeremy L Thompson if (ierr) return ierr; 931ccaff030SJeremy L Thompson 932ccaff030SJeremy L Thompson // Allocate PETSc context 933ccaff030SJeremy L Thompson ierr = PetscCalloc1(1, &user); CHKERRQ(ierr); 934ccaff030SJeremy L Thompson ierr = PetscMalloc1(1, &units); CHKERRQ(ierr); 935ccaff030SJeremy L Thompson 936ccaff030SJeremy L Thompson // Parse command line options 937ccaff030SJeremy L Thompson comm = PETSC_COMM_WORLD; 938ccaff030SJeremy L Thompson ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", 939ccaff030SJeremy L Thompson NULL); CHKERRQ(ierr); 940ccaff030SJeremy L Thompson ierr = PetscOptionsString("-ceed", "CEED resource specifier", 941ccaff030SJeremy L Thompson NULL, ceedresource, ceedresource, 942ccaff030SJeremy L Thompson sizeof(ceedresource), NULL); CHKERRQ(ierr); 943*dc8efd83SLeila Ghaffari ierr = PetscOptionsBool("-test", "Run in test mode", 944*dc8efd83SLeila Ghaffari NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr); 945*dc8efd83SLeila Ghaffari ierr = PetscOptionsScalar("-compare_final_state_atol", 946*dc8efd83SLeila Ghaffari "Test absolute tolerance", 947*dc8efd83SLeila Ghaffari NULL, testtol, &testtol, NULL); CHKERRQ(ierr); 948*dc8efd83SLeila Ghaffari ierr = PetscOptionsString("-compare_final_state_filename", "Test filename", 949*dc8efd83SLeila Ghaffari NULL, filepath, filepath, 950*dc8efd83SLeila Ghaffari sizeof(filepath), NULL); CHKERRQ(ierr); 951ccaff030SJeremy L Thompson problemChoice = NS_DENSITY_CURRENT; 952ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL, 953ccaff030SJeremy L Thompson problemTypes, (PetscEnum)problemChoice, 954ccaff030SJeremy L Thompson (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr); 955ccaff030SJeremy L Thompson problem = &problemOptions[problemChoice]; 9561184866aSLeila Ghaffari ierr = PetscOptionsEnum("-problem_advection_wind", "Wind type in Advection", 957f259b054Svaleriabarra NULL, WindTypes, 958f259b054Svaleriabarra (PetscEnum)(wind_type = ADVECTION_WIND_ROTATION), 9591184866aSLeila Ghaffari (PetscEnum *)&wind_type, NULL); CHKERRQ(ierr); 9601184866aSLeila Ghaffari PetscInt n = problem->dim; 96182c09b01SLeila Ghaffari PetscBool userWind; 962ebb4b9bdSLeila Ghaffari ierr = PetscOptionsRealArray("-problem_advection_wind_translation", 963ebb4b9bdSLeila Ghaffari "Constant wind vector", 96482c09b01SLeila Ghaffari NULL, wind, &n, &userWind); CHKERRQ(ierr); 96582c09b01SLeila Ghaffari if (wind_type == ADVECTION_WIND_ROTATION && userWind) { 966ebb4b9bdSLeila Ghaffari ierr = PetscPrintf(comm, 967ebb4b9bdSLeila Ghaffari "Warning! Use -problem_advection_wind_translation only with -problem_advection_wind translation\n"); 96882c09b01SLeila Ghaffari CHKERRQ(ierr); 9691184866aSLeila Ghaffari } 970ebb4b9bdSLeila Ghaffari if (wind_type == ADVECTION_WIND_TRANSLATION 971ebb4b9bdSLeila Ghaffari && problemChoice == NS_DENSITY_CURRENT) { 97267babd9cSLeila Ghaffari SETERRQ(comm, PETSC_ERR_ARG_INCOMP, 97367babd9cSLeila Ghaffari "-problem_advection_wind translation is not defined for -problem density_current"); 97467babd9cSLeila Ghaffari } 975ccaff030SJeremy L Thompson ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL, 976ccaff030SJeremy L Thompson StabilizationTypes, (PetscEnum)(stab = STAB_NONE), 977ccaff030SJeremy L Thompson (PetscEnum *)&stab, NULL); CHKERRQ(ierr); 978ccaff030SJeremy L Thompson ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation", 979ccaff030SJeremy L Thompson NULL, implicit=PETSC_FALSE, &implicit, NULL); 980ccaff030SJeremy L Thompson CHKERRQ(ierr); 98184d34d69SLeila Ghaffari if (!implicit && stab != STAB_NONE) { 98284d34d69SLeila Ghaffari ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n"); 98384d34d69SLeila Ghaffari CHKERRQ(ierr); 98484d34d69SLeila Ghaffari } 985ccaff030SJeremy L Thompson { 9867573aee6SJed Brown PetscInt len; 9877573aee6SJed Brown PetscBool flg; 988ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray("-bc_wall", 989ccaff030SJeremy L Thompson "Use wall boundary conditions on this list of faces", 990ccaff030SJeremy L Thompson NULL, bc.walls, 991ccaff030SJeremy L Thompson (len = sizeof(bc.walls) / sizeof(bc.walls[0]), 992ccaff030SJeremy L Thompson &len), &flg); CHKERRQ(ierr); 9937573aee6SJed Brown if (flg) { 9947573aee6SJed Brown bc.nwall = len; 9957573aee6SJed Brown // Using a no-slip wall disables automatic slip walls (they must be set explicitly) 9967573aee6SJed Brown bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0; 9977573aee6SJed Brown } 998ccaff030SJeremy L Thompson for (PetscInt j=0; j<3; j++) { 999ccaff030SJeremy L Thompson const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 1000ccaff030SJeremy L Thompson ierr = PetscOptionsIntArray(flags[j], 1001ccaff030SJeremy L Thompson "Use slip boundary conditions on this list of faces", 1002ccaff030SJeremy L Thompson NULL, bc.slips[j], 1003ccaff030SJeremy L Thompson (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]), 1004ccaff030SJeremy L Thompson &len), &flg); 1005ccaff030SJeremy L Thompson CHKERRQ(ierr); 100684d34d69SLeila Ghaffari if (flg) { 100784d34d69SLeila Ghaffari bc.nslip[j] = len; 100884d34d69SLeila Ghaffari bc.userbc = PETSC_TRUE; 100984d34d69SLeila Ghaffari } 1010ccaff030SJeremy L Thompson } 1011ccaff030SJeremy L Thompson } 1012cb3e2689Svaleriabarra ierr = PetscOptionsInt("-viz_refine", 1013cb3e2689Svaleriabarra "Regular refinement levels for visualization", 1014cb3e2689Svaleriabarra NULL, viz_refine, &viz_refine, NULL); 1015ccaff030SJeremy L Thompson CHKERRQ(ierr); 1016ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units", 1017ccaff030SJeremy L Thompson NULL, meter, &meter, NULL); CHKERRQ(ierr); 1018ccaff030SJeremy L Thompson meter = fabs(meter); 1019ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_second","1 second in scaled time units", 1020ccaff030SJeremy L Thompson NULL, second, &second, NULL); CHKERRQ(ierr); 1021ccaff030SJeremy L Thompson second = fabs(second); 1022ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units", 1023ccaff030SJeremy L Thompson NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr); 1024ccaff030SJeremy L Thompson kilogram = fabs(kilogram); 1025ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-units_Kelvin", 1026ccaff030SJeremy L Thompson "1 Kelvin in scaled temperature units", 1027ccaff030SJeremy L Thompson NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr); 1028ccaff030SJeremy L Thompson Kelvin = fabs(Kelvin); 1029ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-theta0", "Reference potential temperature", 1030ccaff030SJeremy L Thompson NULL, theta0, &theta0, NULL); CHKERRQ(ierr); 1031ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature", 1032ccaff030SJeremy L Thompson NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr); 1033ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-P0", "Atmospheric pressure", 1034ccaff030SJeremy L Thompson NULL, P0, &P0, NULL); CHKERRQ(ierr); 103516c0476cSLeila Ghaffari ierr = PetscOptionsScalar("-E_wind", "Total energy of inflow wind", 103616c0476cSLeila Ghaffari NULL, E_wind, &E_wind, NULL); CHKERRQ(ierr); 1037ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency", 1038ccaff030SJeremy L Thompson NULL, N, &N, NULL); CHKERRQ(ierr); 1039ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume", 1040ccaff030SJeremy L Thompson NULL, cv, &cv, NULL); CHKERRQ(ierr); 1041ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure", 1042ccaff030SJeremy L Thompson NULL, cp, &cp, NULL); CHKERRQ(ierr); 1043ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-g", "Gravitational acceleration", 1044ccaff030SJeremy L Thompson NULL, g, &g, NULL); CHKERRQ(ierr); 1045ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lambda", 1046ccaff030SJeremy L Thompson "Stokes hypothesis second viscosity coefficient", 1047ccaff030SJeremy L Thompson NULL, lambda, &lambda, NULL); CHKERRQ(ierr); 1048ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient", 1049ccaff030SJeremy L Thompson NULL, mu, &mu, NULL); CHKERRQ(ierr); 1050ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-k", "Thermal conductivity", 1051ccaff030SJeremy L Thompson NULL, k, &k, NULL); CHKERRQ(ierr); 1052ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-CtauS", 1053ccaff030SJeremy L Thompson "Scale coefficient for tau (nondimensional)", 1054ccaff030SJeremy L Thompson NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr); 105584d34d69SLeila Ghaffari if (stab == STAB_NONE && CtauS != 0) { 105684d34d69SLeila Ghaffari ierr = PetscPrintf(comm, 105784d34d69SLeila Ghaffari "Warning! Use -CtauS only with -stab su or -stab supg\n"); 105884d34d69SLeila Ghaffari CHKERRQ(ierr); 105984d34d69SLeila Ghaffari } 1060ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-strong_form", 1061ccaff030SJeremy L Thompson "Strong (1) or weak/integrated by parts (0) advection residual", 1062ccaff030SJeremy L Thompson NULL, strong_form, &strong_form, NULL); 1063ccaff030SJeremy L Thompson CHKERRQ(ierr); 106484d34d69SLeila Ghaffari if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) { 106584d34d69SLeila Ghaffari ierr = PetscPrintf(comm, 106684d34d69SLeila Ghaffari "Warning! Problem density_current does not support -CtauS or -strong_form\n"); 106784d34d69SLeila Ghaffari CHKERRQ(ierr); 106884d34d69SLeila Ghaffari } 1069ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lx", "Length scale in x direction", 1070ccaff030SJeremy L Thompson NULL, lx, &lx, NULL); CHKERRQ(ierr); 1071ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-ly", "Length scale in y direction", 1072ccaff030SJeremy L Thompson NULL, ly, &ly, NULL); CHKERRQ(ierr); 1073ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-lz", "Length scale in z direction", 1074ccaff030SJeremy L Thompson NULL, lz, &lz, NULL); CHKERRQ(ierr); 1075ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble", 1076ccaff030SJeremy L Thompson NULL, rc, &rc, NULL); CHKERRQ(ierr); 1077ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resx","Target resolution in x", 1078ccaff030SJeremy L Thompson NULL, resx, &resx, NULL); CHKERRQ(ierr); 1079ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resy","Target resolution in y", 1080ccaff030SJeremy L Thompson NULL, resy, &resy, NULL); CHKERRQ(ierr); 1081ccaff030SJeremy L Thompson ierr = PetscOptionsScalar("-resz","Target resolution in z", 1082ccaff030SJeremy L Thompson NULL, resz, &resz, NULL); CHKERRQ(ierr); 108382c09b01SLeila Ghaffari n = problem->dim; 1084ccaff030SJeremy L Thompson center[0] = 0.5 * lx; 1085ccaff030SJeremy L Thompson center[1] = 0.5 * ly; 1086ccaff030SJeremy L Thompson center[2] = 0.5 * lz; 1087ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-center", "Location of bubble center", 1088ccaff030SJeremy L Thompson NULL, center, &n, NULL); CHKERRQ(ierr); 1089ccaff030SJeremy L Thompson n = problem->dim; 1090ccaff030SJeremy L Thompson ierr = PetscOptionsRealArray("-dc_axis", 1091ccaff030SJeremy L Thompson "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric", 1092ccaff030SJeremy L Thompson NULL, dc_axis, &n, NULL); CHKERRQ(ierr); 1093ccaff030SJeremy L Thompson { 1094ccaff030SJeremy L Thompson PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) + 1095ccaff030SJeremy L Thompson PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2])); 1096ccaff030SJeremy L Thompson if (norm > 0) { 1097ccaff030SJeremy L Thompson for (int i=0; i<3; i++) dc_axis[i] /= norm; 1098ccaff030SJeremy L Thompson } 1099ccaff030SJeremy L Thompson } 1100ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-output_freq", 1101ccaff030SJeremy L Thompson "Frequency of output, in number of steps", 1102ccaff030SJeremy L Thompson NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr); 1103ccaff030SJeremy L Thompson ierr = PetscOptionsInt("-continue", "Continue from previous solution", 1104ccaff030SJeremy L Thompson NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr); 110584d34d69SLeila Ghaffari ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements", 110684d34d69SLeila Ghaffari NULL, degree, °ree, NULL); CHKERRQ(ierr); 110784d34d69SLeila Ghaffari ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points", 110884d34d69SLeila Ghaffari NULL, qextra, &qextra, NULL); CHKERRQ(ierr); 110981f92cf0SLeila Ghaffari PetscBool userQextraSur; 1110ebb4b9bdSLeila Ghaffari ierr = PetscOptionsInt("-qextra_boundary", 1111ebb4b9bdSLeila Ghaffari "Number of extra quadrature points on in/outflow faces", 1112f259b054Svaleriabarra NULL, qextraSur, &qextraSur, &userQextraSur); 1113f259b054Svaleriabarra CHKERRQ(ierr); 1114ebb4b9bdSLeila Ghaffari if ((wind_type == ADVECTION_WIND_ROTATION 1115ebb4b9bdSLeila Ghaffari || problemChoice == NS_DENSITY_CURRENT) && userQextraSur) { 1116ebb4b9bdSLeila Ghaffari ierr = PetscPrintf(comm, 1117ebb4b9bdSLeila Ghaffari "Warning! Use -qextra_boundary only with -problem_advection_wind translation\n"); 111881f92cf0SLeila Ghaffari CHKERRQ(ierr); 111981f92cf0SLeila Ghaffari } 1120d99129b9SLeila Ghaffari ierr = PetscStrncpy(user->outputdir, ".", 2); CHKERRQ(ierr); 1121d99129b9SLeila Ghaffari ierr = PetscOptionsString("-output_dir", "Output directory", 1122d99129b9SLeila Ghaffari NULL, user->outputdir, user->outputdir, 1123d99129b9SLeila Ghaffari sizeof(user->outputdir), NULL); CHKERRQ(ierr); 112484d34d69SLeila Ghaffari memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST; 112584d34d69SLeila Ghaffari ierr = PetscOptionsEnum("-memtype", 112684d34d69SLeila Ghaffari "CEED MemType requested", NULL, 112784d34d69SLeila Ghaffari memTypes, (PetscEnum)memtyperequested, 112884d34d69SLeila Ghaffari (PetscEnum *)&memtyperequested, &setmemtyperequest); 112984d34d69SLeila Ghaffari CHKERRQ(ierr); 1130ccaff030SJeremy L Thompson ierr = PetscOptionsEnd(); CHKERRQ(ierr); 1131ccaff030SJeremy L Thompson 1132ccaff030SJeremy L Thompson // Define derived units 1133ccaff030SJeremy L Thompson Pascal = kilogram / (meter * PetscSqr(second)); 1134ccaff030SJeremy L Thompson JperkgK = PetscSqr(meter) / (PetscSqr(second) * Kelvin); 1135ccaff030SJeremy L Thompson mpersquareds = meter / PetscSqr(second); 1136ccaff030SJeremy L Thompson WpermK = kilogram * meter / (pow(second,3) * Kelvin); 1137ccaff030SJeremy L Thompson kgpercubicm = kilogram / pow(meter,3); 1138ccaff030SJeremy L Thompson kgpersquaredms = kilogram / (PetscSqr(meter) * second); 1139ccaff030SJeremy L Thompson Joulepercubicm = kilogram / (meter * PetscSqr(second)); 114016c0476cSLeila Ghaffari Joule = kilogram * PetscSqr(meter) / PetscSqr(second); 1141ccaff030SJeremy L Thompson 1142ccaff030SJeremy L Thompson // Scale variables to desired units 1143ccaff030SJeremy L Thompson theta0 *= Kelvin; 1144ccaff030SJeremy L Thompson thetaC *= Kelvin; 1145ccaff030SJeremy L Thompson P0 *= Pascal; 114616c0476cSLeila Ghaffari E_wind *= Joule; 1147ccaff030SJeremy L Thompson N *= (1./second); 1148ccaff030SJeremy L Thompson cv *= JperkgK; 1149ccaff030SJeremy L Thompson cp *= JperkgK; 1150ccaff030SJeremy L Thompson Rd = cp - cv; 1151ccaff030SJeremy L Thompson g *= mpersquareds; 1152ccaff030SJeremy L Thompson mu *= Pascal * second; 1153ccaff030SJeremy L Thompson k *= WpermK; 1154ccaff030SJeremy L Thompson lx = fabs(lx) * meter; 1155ccaff030SJeremy L Thompson ly = fabs(ly) * meter; 1156ccaff030SJeremy L Thompson lz = fabs(lz) * meter; 1157ccaff030SJeremy L Thompson rc = fabs(rc) * meter; 1158ccaff030SJeremy L Thompson resx = fabs(resx) * meter; 1159ccaff030SJeremy L Thompson resy = fabs(resy) * meter; 1160ccaff030SJeremy L Thompson resz = fabs(resz) * meter; 1161ccaff030SJeremy L Thompson for (int i=0; i<3; i++) center[i] *= meter; 1162ccaff030SJeremy L Thompson 1163ccaff030SJeremy L Thompson const CeedInt dim = problem->dim, ncompx = problem->dim, 1164cfa64770SLeila Ghaffari qdatasizeVol = problem->qdatasizeVol; 1165ccaff030SJeremy L Thompson // Set up the libCEED context 1166777ff853SJeremy L Thompson struct SetupContext_ ctxSetupData = { 1167ccaff030SJeremy L Thompson .theta0 = theta0, 1168ccaff030SJeremy L Thompson .thetaC = thetaC, 1169ccaff030SJeremy L Thompson .P0 = P0, 1170ccaff030SJeremy L Thompson .N = N, 1171ccaff030SJeremy L Thompson .cv = cv, 1172ccaff030SJeremy L Thompson .cp = cp, 1173ccaff030SJeremy L Thompson .Rd = Rd, 1174ccaff030SJeremy L Thompson .g = g, 1175ccaff030SJeremy L Thompson .rc = rc, 1176ccaff030SJeremy L Thompson .lx = lx, 1177ccaff030SJeremy L Thompson .ly = ly, 1178ccaff030SJeremy L Thompson .lz = lz, 1179ccaff030SJeremy L Thompson .center[0] = center[0], 1180ccaff030SJeremy L Thompson .center[1] = center[1], 1181ccaff030SJeremy L Thompson .center[2] = center[2], 1182ccaff030SJeremy L Thompson .dc_axis[0] = dc_axis[0], 1183ccaff030SJeremy L Thompson .dc_axis[1] = dc_axis[1], 1184ccaff030SJeremy L Thompson .dc_axis[2] = dc_axis[2], 11851184866aSLeila Ghaffari .wind[0] = wind[0], 11861184866aSLeila Ghaffari .wind[1] = wind[1], 11871184866aSLeila Ghaffari .wind[2] = wind[2], 1188ccaff030SJeremy L Thompson .time = 0, 11891184866aSLeila Ghaffari .wind_type = wind_type, 1190ccaff030SJeremy L Thompson }; 1191ccaff030SJeremy L Thompson 119284d34d69SLeila Ghaffari // Create the mesh 1193ccaff030SJeremy L Thompson { 1194ccaff030SJeremy L Thompson const PetscReal scale[3] = {lx, ly, lz}; 1195ccaff030SJeremy L Thompson ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale, 119684d34d69SLeila Ghaffari NULL, PETSC_TRUE, &dm); 1197ccaff030SJeremy L Thompson CHKERRQ(ierr); 1198ccaff030SJeremy L Thompson } 119984d34d69SLeila Ghaffari 120084d34d69SLeila Ghaffari // Distribute the mesh over processes 120184d34d69SLeila Ghaffari { 1202ccaff030SJeremy L Thompson DM dmDist = NULL; 1203ccaff030SJeremy L Thompson PetscPartitioner part; 1204ccaff030SJeremy L Thompson 1205ccaff030SJeremy L Thompson ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr); 1206ccaff030SJeremy L Thompson ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr); 1207ccaff030SJeremy L Thompson ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr); 1208ccaff030SJeremy L Thompson if (dmDist) { 1209ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 1210ccaff030SJeremy L Thompson dm = dmDist; 1211ccaff030SJeremy L Thompson } 1212ccaff030SJeremy L Thompson } 1213ccaff030SJeremy L Thompson ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr); 1214ccaff030SJeremy L Thompson 121584d34d69SLeila Ghaffari // Setup DM 1216ccaff030SJeremy L Thompson ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr); 1217ccaff030SJeremy L Thompson ierr = DMSetFromOptions(dm); CHKERRQ(ierr); 1218777ff853SJeremy L Thompson ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetupData); CHKERRQ(ierr); 121984d34d69SLeila Ghaffari 122084d34d69SLeila Ghaffari // Refine DM for high-order viz 1221ccaff030SJeremy L Thompson dmviz = NULL; 1222ccaff030SJeremy L Thompson interpviz = NULL; 1223ccaff030SJeremy L Thompson if (viz_refine) { 1224ff6701fcSJed Brown DM dmhierarchy[viz_refine+1]; 1225ff6701fcSJed Brown 1226ccaff030SJeremy L Thompson ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr); 1227ff6701fcSJed Brown dmhierarchy[0] = dm; 122884d34d69SLeila Ghaffari for (PetscInt i = 0, d = degree; i < viz_refine; i++) { 1229ff6701fcSJed Brown Mat interp_next; 1230ff6701fcSJed Brown 1231ff6701fcSJed Brown ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]); 1232ccaff030SJeremy L Thompson CHKERRQ(ierr); 1233bc6a41f7SJed Brown ierr = DMClearDS(dmhierarchy[i+1]); CHKERRQ(ierr); 1234bc6a41f7SJed Brown ierr = DMClearFields(dmhierarchy[i+1]); CHKERRQ(ierr); 1235ff6701fcSJed Brown ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr); 1236ff6701fcSJed Brown d = (d + 1) / 2; 1237ff6701fcSJed Brown if (i + 1 == viz_refine) d = 1; 1238777ff853SJeremy L Thompson ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetupData); 1239f259b054Svaleriabarra CHKERRQ(ierr); 1240ff6701fcSJed Brown ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1], 1241ff6701fcSJed Brown &interp_next, NULL); CHKERRQ(ierr); 1242ff6701fcSJed Brown if (!i) interpviz = interp_next; 1243ff6701fcSJed Brown else { 1244ff6701fcSJed Brown Mat C; 1245ff6701fcSJed Brown ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX, 1246ff6701fcSJed Brown PETSC_DECIDE, &C); CHKERRQ(ierr); 1247ff6701fcSJed Brown ierr = MatDestroy(&interp_next); CHKERRQ(ierr); 1248ff6701fcSJed Brown ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1249ff6701fcSJed Brown interpviz = C; 1250ff6701fcSJed Brown } 1251ff6701fcSJed Brown } 1252cb3e2689Svaleriabarra for (PetscInt i=1; i<viz_refine; i++) { 1253ff6701fcSJed Brown ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr); 1254cb3e2689Svaleriabarra } 1255ff6701fcSJed Brown dmviz = dmhierarchy[viz_refine]; 1256ccaff030SJeremy L Thompson } 1257ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr); 1258ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr); 1259ccaff030SJeremy L Thompson ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr); 1260ccaff030SJeremy L Thompson lnodes /= ncompq; 1261ccaff030SJeremy L Thompson 126284d34d69SLeila Ghaffari // Initialize CEED 126384d34d69SLeila Ghaffari CeedInit(ceedresource, &ceed); 126484d34d69SLeila Ghaffari // Set memtype 126584d34d69SLeila Ghaffari CeedMemType memtypebackend; 126684d34d69SLeila Ghaffari CeedGetPreferredMemType(ceed, &memtypebackend); 126784d34d69SLeila Ghaffari // Check memtype compatibility 126884d34d69SLeila Ghaffari if (!setmemtyperequest) 126984d34d69SLeila Ghaffari memtyperequested = memtypebackend; 127084d34d69SLeila Ghaffari else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE) 127184d34d69SLeila Ghaffari SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS, 127284d34d69SLeila Ghaffari "PETSc was not built with CUDA. " 127384d34d69SLeila Ghaffari "Requested MemType CEED_MEM_DEVICE is not supported.", NULL); 127484d34d69SLeila Ghaffari 127584d34d69SLeila Ghaffari // Set number of 1D nodes and quadrature points 127684d34d69SLeila Ghaffari numP = degree + 1; 127784d34d69SLeila Ghaffari numQ = numP + qextra; 127884d34d69SLeila Ghaffari 127984d34d69SLeila Ghaffari // Print summary 1280*dc8efd83SLeila Ghaffari if (!test) { 1281ccaff030SJeremy L Thompson CeedInt gdofs, odofs; 1282ccaff030SJeremy L Thompson int comm_size; 1283ccaff030SJeremy L Thompson char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE"; 1284ccaff030SJeremy L Thompson ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr); 1285ccaff030SJeremy L Thompson ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr); 128684d34d69SLeila Ghaffari gnodes = gdofs/ncompq; 1287ccaff030SJeremy L Thompson ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr); 1288ccaff030SJeremy L Thompson ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, 1289ccaff030SJeremy L Thompson sizeof(box_faces_str), NULL); CHKERRQ(ierr); 129084d34d69SLeila Ghaffari const char *usedresource; 129184d34d69SLeila Ghaffari CeedGetResource(ceed, &usedresource); 1292ccaff030SJeremy L Thompson 129384d34d69SLeila Ghaffari ierr = PetscPrintf(comm, 129484d34d69SLeila Ghaffari "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 129584d34d69SLeila Ghaffari " rank(s) : %d\n" 129684d34d69SLeila Ghaffari " Problem:\n" 129784d34d69SLeila Ghaffari " Problem Name : %s\n" 129884d34d69SLeila Ghaffari " Stabilization : %s\n" 129984d34d69SLeila Ghaffari " PETSc:\n" 130084d34d69SLeila Ghaffari " Box Faces : %s\n" 130184d34d69SLeila Ghaffari " libCEED:\n" 130284d34d69SLeila Ghaffari " libCEED Backend : %s\n" 130384d34d69SLeila Ghaffari " libCEED Backend MemType : %s\n" 130484d34d69SLeila Ghaffari " libCEED User Requested MemType : %s\n" 130584d34d69SLeila Ghaffari " Mesh:\n" 130684d34d69SLeila Ghaffari " Number of 1D Basis Nodes (P) : %d\n" 130784d34d69SLeila Ghaffari " Number of 1D Quadrature Points (Q) : %d\n" 130884d34d69SLeila Ghaffari " Global DoFs : %D\n" 130984d34d69SLeila Ghaffari " Owned DoFs : %D\n" 131084d34d69SLeila Ghaffari " DoFs per node : %D\n" 131184d34d69SLeila Ghaffari " Global nodes : %D\n" 131284d34d69SLeila Ghaffari " Owned nodes : %D\n", 131384d34d69SLeila Ghaffari comm_size, problemTypes[problemChoice], 131484d34d69SLeila Ghaffari StabilizationTypes[stab], box_faces_str, usedresource, 131584d34d69SLeila Ghaffari CeedMemTypes[memtypebackend], 131684d34d69SLeila Ghaffari (setmemtyperequest) ? 131784d34d69SLeila Ghaffari CeedMemTypes[memtyperequested] : "none", 131884d34d69SLeila Ghaffari numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes); 131984d34d69SLeila Ghaffari CHKERRQ(ierr); 13200c6c0b13SLeila Ghaffari } 13210c6c0b13SLeila Ghaffari 1322ccaff030SJeremy L Thompson // Set up global mass vector 1323ccaff030SJeremy L Thompson ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr); 1324ccaff030SJeremy L Thompson 132584d34d69SLeila Ghaffari // Set up libCEED 1326ccaff030SJeremy L Thompson // CEED Bases 1327ccaff030SJeremy L Thompson CeedInit(ceedresource, &ceed); 132884d34d69SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS, 132984d34d69SLeila Ghaffari &basisq); 133084d34d69SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS, 133184d34d69SLeila Ghaffari &basisx); 133284d34d69SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP, 133384d34d69SLeila Ghaffari CEED_GAUSS_LOBATTO, &basisxc); 1334ccaff030SJeremy L Thompson ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr); 1335ccaff030SJeremy L Thompson ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL); 1336ccaff030SJeremy L Thompson CHKERRQ(ierr); 1337ccaff030SJeremy L Thompson 1338ccaff030SJeremy L Thompson // CEED Restrictions 13391e150236SLeila Ghaffari ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ, 134084d34d69SLeila Ghaffari qdatasizeVol, &restrictq, &restrictx, 134184d34d69SLeila Ghaffari &restrictqdi); CHKERRQ(ierr); 1342ccaff030SJeremy L Thompson 1343ccaff030SJeremy L Thompson ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr); 1344ccaff030SJeremy L Thompson ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr); 1345ccaff030SJeremy L Thompson 1346ccaff030SJeremy L Thompson // Create the CEED vectors that will be needed in setup 1347bd910870SLeila Ghaffari CeedInt NqptsVol; 134884d34d69SLeila Ghaffari CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol); 134984d34d69SLeila Ghaffari CeedElemRestrictionGetNumElements(restrictq, &localNelemVol); 13508b982baeSLeila Ghaffari CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata); 135184d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL); 1352ccaff030SJeremy L Thompson 1353ccaff030SJeremy L Thompson // Create the Q-function that builds the quadrature data for the NS operator 1354ea6e0f84SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc, 1355ea6e0f84SLeila Ghaffari &qf_setupVol); 1356ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD); 1357ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT); 13588b982baeSLeila Ghaffari CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1359ccaff030SJeremy L Thompson 1360ccaff030SJeremy L Thompson // Create the Q-function that sets the ICs of the operator 1361ccaff030SJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics); 1362ccaff030SJeremy L Thompson CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP); 1363ccaff030SJeremy L Thompson CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE); 1364ccaff030SJeremy L Thompson 1365ea6e0f84SLeila Ghaffari qf_rhsVol = NULL; 1366ea6e0f84SLeila Ghaffari if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator 1367ea6e0f84SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs, 1368ea6e0f84SLeila Ghaffari problem->applyVol_rhs_loc, &qf_rhsVol); 1369ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP); 1370ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD); 13718b982baeSLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1372ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP); 1373ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP); 1374ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD); 1375ccaff030SJeremy L Thompson } 1376ccaff030SJeremy L Thompson 1377ea6e0f84SLeila Ghaffari qf_ifunctionVol = NULL; 1378ea6e0f84SLeila Ghaffari if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction 1379ea6e0f84SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction, 1380ea6e0f84SLeila Ghaffari problem->applyVol_ifunction_loc, &qf_ifunctionVol); 1381ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP); 1382ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD); 1383ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP); 13848b982baeSLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE); 1385ea6e0f84SLeila Ghaffari CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP); 1386ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP); 1387ea6e0f84SLeila Ghaffari CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD); 1388ccaff030SJeremy L Thompson } 1389ccaff030SJeremy L Thompson 1390ccaff030SJeremy L Thompson // Create the operator that builds the quadrature data for the NS operator 1391ea6e0f84SLeila Ghaffari CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol); 139284d34d69SLeila Ghaffari CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE); 1393ea6e0f84SLeila Ghaffari CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE, 139484d34d69SLeila Ghaffari basisx, CEED_VECTOR_NONE); 139584d34d69SLeila Ghaffari CeedOperatorSetField(op_setupVol, "qdata", restrictqdi, 1396ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1397ccaff030SJeremy L Thompson 1398ccaff030SJeremy L Thompson // Create the operator that sets the ICs 1399ccaff030SJeremy L Thompson CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics); 140084d34d69SLeila Ghaffari CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE); 140184d34d69SLeila Ghaffari CeedOperatorSetField(op_ics, "q0", restrictq, 1402ccaff030SJeremy L Thompson CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 1403ccaff030SJeremy L Thompson 140484d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL); 140584d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL); 140684d34d69SLeila Ghaffari CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL); 1407ccaff030SJeremy L Thompson 1408ea6e0f84SLeila Ghaffari if (qf_rhsVol) { // Create the RHS physics operator 1409ccaff030SJeremy L Thompson CeedOperator op; 1410ea6e0f84SLeila Ghaffari CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op); 141184d34d69SLeila Ghaffari CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 141284d34d69SLeila Ghaffari CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 141384d34d69SLeila Ghaffari CeedOperatorSetField(op, "qdata", restrictqdi, 14148b982baeSLeila Ghaffari CEED_BASIS_COLLOCATED, qdata); 141584d34d69SLeila Ghaffari CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 141684d34d69SLeila Ghaffari CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 141784d34d69SLeila Ghaffari CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1418d3630711SJed Brown user->op_rhs_vol = op; 1419ccaff030SJeremy L Thompson } 1420ccaff030SJeremy L Thompson 1421ea6e0f84SLeila Ghaffari if (qf_ifunctionVol) { // Create the IFunction operator 1422ccaff030SJeremy L Thompson CeedOperator op; 1423ea6e0f84SLeila Ghaffari CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op); 142484d34d69SLeila Ghaffari CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE); 142584d34d69SLeila Ghaffari CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE); 142684d34d69SLeila Ghaffari CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed); 142784d34d69SLeila Ghaffari CeedOperatorSetField(op, "qdata", restrictqdi, 14288b982baeSLeila Ghaffari CEED_BASIS_COLLOCATED, qdata); 142984d34d69SLeila Ghaffari CeedOperatorSetField(op, "x", restrictx, basisx, xcorners); 143084d34d69SLeila Ghaffari CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE); 143184d34d69SLeila Ghaffari CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE); 1432d3630711SJed Brown user->op_ifunction_vol = op; 1433ccaff030SJeremy L Thompson } 1434ccaff030SJeremy L Thompson 14356a0edaf9SLeila Ghaffari // Set up CEED for the boundaries 14366a0edaf9SLeila Ghaffari CeedInt height = 1; 14376a0edaf9SLeila Ghaffari CeedInt dimSur = dim - height; 14381e150236SLeila Ghaffari CeedInt numP_Sur = degree + 1; 14391e150236SLeila Ghaffari CeedInt numQ_Sur = numP_Sur + qextraSur; 1440cfa64770SLeila Ghaffari const CeedInt qdatasizeSur = problem->qdatasizeSur; 1441cfa64770SLeila Ghaffari CeedBasis basisxSur, basisxcSur, basisqSur; 1442cfa64770SLeila Ghaffari CeedInt NqptsSur; 14437ed8b9c7SLeila Ghaffari CeedQFunction qf_setupSur, qf_applySur; 1444cfa64770SLeila Ghaffari 1445cfa64770SLeila Ghaffari // CEED bases for the boundaries 1446ebb4b9bdSLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, 1447ebb4b9bdSLeila Ghaffari CEED_GAUSS, 14486a0edaf9SLeila Ghaffari &basisqSur); 14496a0edaf9SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS, 14506a0edaf9SLeila Ghaffari &basisxSur); 14516a0edaf9SLeila Ghaffari CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur, 14526a0edaf9SLeila Ghaffari CEED_GAUSS_LOBATTO, &basisxcSur); 14536a0edaf9SLeila Ghaffari CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur); 14546a0edaf9SLeila Ghaffari 1455cfa64770SLeila Ghaffari // Create the Q-function that builds the quadrature data for the Surface operator 14566a0edaf9SLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc, 14576a0edaf9SLeila Ghaffari &qf_setupSur); 14586a0edaf9SLeila Ghaffari CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD); 14596a0edaf9SLeila Ghaffari CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT); 14606a0edaf9SLeila Ghaffari CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 14616a0edaf9SLeila Ghaffari 14627659d40cSLeila Ghaffari // Creat Q-Function for Boundaries 14637ed8b9c7SLeila Ghaffari qf_applySur = NULL; 14647659d40cSLeila Ghaffari if (problem->applySur) { 14657659d40cSLeila Ghaffari CeedQFunctionCreateInterior(ceed, 1, problem->applySur, 14667ed8b9c7SLeila Ghaffari problem->applySur_loc, &qf_applySur); 14677ed8b9c7SLeila Ghaffari CeedQFunctionAddInput(qf_applySur, "q", ncompq, CEED_EVAL_INTERP); 14687ed8b9c7SLeila Ghaffari CeedQFunctionAddInput(qf_applySur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE); 14697ed8b9c7SLeila Ghaffari CeedQFunctionAddInput(qf_applySur, "x", ncompx, CEED_EVAL_INTERP); 14707ed8b9c7SLeila Ghaffari CeedQFunctionAddOutput(qf_applySur, "v", ncompq, CEED_EVAL_INTERP); 14719fe13df9SLeila Ghaffari } 14729fe13df9SLeila Ghaffari 14739fe13df9SLeila Ghaffari // Create CEED Operator for the whole domain 14749fe13df9SLeila Ghaffari if (!implicit) 1475ebb4b9bdSLeila Ghaffari ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_rhs_vol, 1476ebb4b9bdSLeila Ghaffari qf_applySur, qf_setupSur, 1477f259b054Svaleriabarra height, numP_Sur, numQ_Sur, qdatasizeSur, 1478f259b054Svaleriabarra NqptsSur, basisxSur, basisqSur, 1479f259b054Svaleriabarra &user->op_rhs); CHKERRQ(ierr); 14809fe13df9SLeila Ghaffari if (implicit) 1481f259b054Svaleriabarra ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, 1482f259b054Svaleriabarra user->op_ifunction_vol, 1483ebb4b9bdSLeila Ghaffari qf_applySur, qf_setupSur, 1484f259b054Svaleriabarra height, numP_Sur, numQ_Sur, qdatasizeSur, 1485f259b054Svaleriabarra NqptsSur, basisxSur, basisqSur, 1486f259b054Svaleriabarra &user->op_ifunction); CHKERRQ(ierr); 14877659d40cSLeila Ghaffari // Set up contex for QFunctions 1488777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxSetup); 1489777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxSetup, CEED_MEM_HOST, CEED_USE_POINTER, 1490777ff853SJeremy L Thompson sizeof ctxSetupData, &ctxSetupData); 1491777ff853SJeremy L Thompson CeedQFunctionSetContext(qf_ics, ctxSetup); 1492777ff853SJeremy L Thompson 1493777ff853SJeremy L Thompson CeedScalar ctxNSData[8] = {lambda, mu, k, cv, cp, g, Rd}; 1494777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxNS); 1495777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxNS, CEED_MEM_HOST, CEED_USE_POINTER, 1496777ff853SJeremy L Thompson sizeof ctxNSData, &ctxNSData); 1497777ff853SJeremy L Thompson 1498777ff853SJeremy L Thompson struct Advection2dContext_ ctxAdvection2dData = { 1499ccaff030SJeremy L Thompson .CtauS = CtauS, 1500ccaff030SJeremy L Thompson .strong_form = strong_form, 1501ccaff030SJeremy L Thompson .stabilization = stab, 1502ccaff030SJeremy L Thompson }; 1503777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxAdvection2d); 1504777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxAdvection2d, CEED_MEM_HOST, CEED_USE_POINTER, 1505777ff853SJeremy L Thompson sizeof ctxAdvection2dData, &ctxAdvection2dData); 1506777ff853SJeremy L Thompson 1507777ff853SJeremy L Thompson struct SurfaceContext_ ctxSurfaceData = { 150816c0476cSLeila Ghaffari .E_wind = E_wind, 15097659d40cSLeila Ghaffari .strong_form = strong_form, 15107659d40cSLeila Ghaffari .implicit = implicit, 15117659d40cSLeila Ghaffari }; 1512777ff853SJeremy L Thompson CeedQFunctionContextCreate(ceed, &ctxSurface); 1513777ff853SJeremy L Thompson CeedQFunctionContextSetData(ctxSurface, CEED_MEM_HOST, CEED_USE_POINTER, 1514777ff853SJeremy L Thompson sizeof ctxSurfaceData, &ctxSurfaceData); 1515777ff853SJeremy L Thompson 1516ccaff030SJeremy L Thompson switch (problemChoice) { 1517ccaff030SJeremy L Thompson case NS_DENSITY_CURRENT: 1518777ff853SJeremy L Thompson if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, ctxNS); 1519777ff853SJeremy L Thompson if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, ctxNS); 1520ccaff030SJeremy L Thompson break; 1521ccaff030SJeremy L Thompson case NS_ADVECTION: 1522ccaff030SJeremy L Thompson case NS_ADVECTION2D: 1523777ff853SJeremy L Thompson if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, ctxAdvection2d); 1524777ff853SJeremy L Thompson if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, ctxAdvection2d); 1525777ff853SJeremy L Thompson if (qf_applySur) CeedQFunctionSetContext(qf_applySur, ctxSurface); 1526ccaff030SJeremy L Thompson } 1527ccaff030SJeremy L Thompson 1528ccaff030SJeremy L Thompson // Set up PETSc context 1529ccaff030SJeremy L Thompson // Set up units structure 1530ccaff030SJeremy L Thompson units->meter = meter; 1531ccaff030SJeremy L Thompson units->kilogram = kilogram; 1532ccaff030SJeremy L Thompson units->second = second; 1533ccaff030SJeremy L Thompson units->Kelvin = Kelvin; 1534ccaff030SJeremy L Thompson units->Pascal = Pascal; 1535ccaff030SJeremy L Thompson units->JperkgK = JperkgK; 1536ccaff030SJeremy L Thompson units->mpersquareds = mpersquareds; 1537ccaff030SJeremy L Thompson units->WpermK = WpermK; 1538ccaff030SJeremy L Thompson units->kgpercubicm = kgpercubicm; 1539ccaff030SJeremy L Thompson units->kgpersquaredms = kgpersquaredms; 1540ccaff030SJeremy L Thompson units->Joulepercubicm = Joulepercubicm; 154116c0476cSLeila Ghaffari units->Joule = Joule; 1542ccaff030SJeremy L Thompson 1543ccaff030SJeremy L Thompson // Set up user structure 1544ccaff030SJeremy L Thompson user->comm = comm; 1545ccaff030SJeremy L Thompson user->outputfreq = outputfreq; 1546ccaff030SJeremy L Thompson user->contsteps = contsteps; 1547ccaff030SJeremy L Thompson user->units = units; 1548ccaff030SJeremy L Thompson user->dm = dm; 1549ccaff030SJeremy L Thompson user->dmviz = dmviz; 1550ccaff030SJeremy L Thompson user->interpviz = interpviz; 1551ccaff030SJeremy L Thompson user->ceed = ceed; 1552ccaff030SJeremy L Thompson 15538b982baeSLeila Ghaffari // Calculate qdata and ICs 1554ccaff030SJeremy L Thompson // Set up state global and local vectors 1555ccaff030SJeremy L Thompson ierr = VecZeroEntries(Q); CHKERRQ(ierr); 1556ccaff030SJeremy L Thompson 1557cfa64770SLeila Ghaffari ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr); 1558ccaff030SJeremy L Thompson 1559ccaff030SJeremy L Thompson // Apply Setup Ceed Operators 1560ccaff030SJeremy L Thompson ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr); 15618b982baeSLeila Ghaffari CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE); 156284d34d69SLeila Ghaffari ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata, 1563ccaff030SJeremy L Thompson user->M); CHKERRQ(ierr); 1564ccaff030SJeremy L Thompson 156584d34d69SLeila Ghaffari ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq, 1566777ff853SJeremy L Thompson ctxSetup, 0.0); CHKERRQ(ierr); 1567ccaff030SJeremy L Thompson if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues() 1568ccaff030SJeremy L Thompson // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow. If we 1569ccaff030SJeremy L Thompson // disable this, we should still get the same results due to the problem->bc function, but with potentially much 1570ccaff030SJeremy L Thompson // slower execution. 1571ccaff030SJeremy L Thompson Vec Qbc; 1572ccaff030SJeremy L Thompson ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1573ccaff030SJeremy L Thompson ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr); 1574ccaff030SJeremy L Thompson ierr = VecZeroEntries(Qloc); CHKERRQ(ierr); 1575ccaff030SJeremy L Thompson ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr); 1576ccaff030SJeremy L Thompson ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr); 1577ccaff030SJeremy L Thompson ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr); 1578ccaff030SJeremy L Thompson ierr = PetscObjectComposeFunction((PetscObject)dm, 157984d34d69SLeila Ghaffari "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS); 158084d34d69SLeila Ghaffari CHKERRQ(ierr); 1581ccaff030SJeremy L Thompson } 1582ccaff030SJeremy L Thompson 1583ccaff030SJeremy L Thompson MPI_Comm_rank(comm, &rank); 1584d99129b9SLeila Ghaffari if (!rank) {ierr = PetscMkdir(user->outputdir); CHKERRQ(ierr);} 1585ccaff030SJeremy L Thompson // Gather initial Q values 1586ccaff030SJeremy L Thompson // In case of continuation of simulation, set up initial values from binary file 1587ccaff030SJeremy L Thompson if (contsteps) { // continue from existent solution 1588ccaff030SJeremy L Thompson PetscViewer viewer; 1589ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1590ccaff030SJeremy L Thompson // Read input 1591ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin", 1592d99129b9SLeila Ghaffari user->outputdir); 1593ccaff030SJeremy L Thompson CHKERRQ(ierr); 1594ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1595ccaff030SJeremy L Thompson CHKERRQ(ierr); 1596ccaff030SJeremy L Thompson ierr = VecLoad(Q, viewer); CHKERRQ(ierr); 1597ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1598ccaff030SJeremy L Thompson } 1599ccaff030SJeremy L Thompson ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr); 1600ccaff030SJeremy L Thompson 1601ccaff030SJeremy L Thompson // Create and setup TS 1602ccaff030SJeremy L Thompson ierr = TSCreate(comm, &ts); CHKERRQ(ierr); 1603ccaff030SJeremy L Thompson ierr = TSSetDM(ts, dm); CHKERRQ(ierr); 1604ccaff030SJeremy L Thompson if (implicit) { 1605ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr); 1606ccaff030SJeremy L Thompson if (user->op_ifunction) { 1607ccaff030SJeremy L Thompson ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr); 1608ccaff030SJeremy L Thompson } else { // Implicit integrators can fall back to using an RHSFunction 1609ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1610ccaff030SJeremy L Thompson } 1611ccaff030SJeremy L Thompson } else { 1612ccaff030SJeremy L Thompson if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL, 1613ccaff030SJeremy L Thompson "Problem does not provide RHSFunction"); 1614ccaff030SJeremy L Thompson ierr = TSSetType(ts, TSRK); CHKERRQ(ierr); 1615ccaff030SJeremy L Thompson ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr); 1616ccaff030SJeremy L Thompson ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr); 1617ccaff030SJeremy L Thompson } 1618ccaff030SJeremy L Thompson ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr); 1619ccaff030SJeremy L Thompson ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr); 1620ccaff030SJeremy L Thompson ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr); 1621*dc8efd83SLeila Ghaffari if (test) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);} 1622ccaff030SJeremy L Thompson ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr); 1623ccaff030SJeremy L Thompson ierr = TSAdaptSetStepLimits(adapt, 1624ccaff030SJeremy L Thompson 1.e-12 * units->second, 1625ccaff030SJeremy L Thompson 1.e2 * units->second); CHKERRQ(ierr); 1626ccaff030SJeremy L Thompson ierr = TSSetFromOptions(ts); CHKERRQ(ierr); 1627ccaff030SJeremy L Thompson if (!contsteps) { // print initial condition 1628*dc8efd83SLeila Ghaffari if (!test) { 1629ccaff030SJeremy L Thompson ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr); 1630ccaff030SJeremy L Thompson } 1631ccaff030SJeremy L Thompson } else { // continue from time of last output 1632ccaff030SJeremy L Thompson PetscReal time; 1633ccaff030SJeremy L Thompson PetscInt count; 1634ccaff030SJeremy L Thompson PetscViewer viewer; 1635ccaff030SJeremy L Thompson char filepath[PETSC_MAX_PATH_LEN]; 1636ccaff030SJeremy L Thompson ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin", 1637d99129b9SLeila Ghaffari user->outputdir); CHKERRQ(ierr); 1638ccaff030SJeremy L Thompson ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 1639ccaff030SJeremy L Thompson CHKERRQ(ierr); 1640ccaff030SJeremy L Thompson ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL); 1641ccaff030SJeremy L Thompson CHKERRQ(ierr); 1642ccaff030SJeremy L Thompson ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 1643ccaff030SJeremy L Thompson ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr); 1644ccaff030SJeremy L Thompson } 1645*dc8efd83SLeila Ghaffari if (!test) { 1646ccaff030SJeremy L Thompson ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr); 1647ccaff030SJeremy L Thompson } 1648ccaff030SJeremy L Thompson 1649ccaff030SJeremy L Thompson // Solve 1650ccaff030SJeremy L Thompson start = MPI_Wtime(); 1651ccaff030SJeremy L Thompson ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr); 1652ccaff030SJeremy L Thompson ierr = TSSolve(ts, Q); CHKERRQ(ierr); 1653ccaff030SJeremy L Thompson cpu_time_used = MPI_Wtime() - start; 1654ccaff030SJeremy L Thompson ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr); 1655ccaff030SJeremy L Thompson ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN, 1656ccaff030SJeremy L Thompson comm); CHKERRQ(ierr); 1657*dc8efd83SLeila Ghaffari if (!test) { 1658ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 165984d34d69SLeila Ghaffari "Time taken for solution (sec): %g\n", 1660ccaff030SJeremy L Thompson (double)cpu_time_used); CHKERRQ(ierr); 1661ccaff030SJeremy L Thompson } 1662ccaff030SJeremy L Thompson 1663ccaff030SJeremy L Thompson // Get error 1664*dc8efd83SLeila Ghaffari if (problem->non_zero_time && !test) { 1665ccaff030SJeremy L Thompson Vec Qexact, Qexactloc; 1666ccaff030SJeremy L Thompson PetscReal norm; 1667ccaff030SJeremy L Thompson ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr); 1668ccaff030SJeremy L Thompson ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 1669ccaff030SJeremy L Thompson ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr); 1670ccaff030SJeremy L Thompson 167184d34d69SLeila Ghaffari ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact, 1672777ff853SJeremy L Thompson restrictq, ctxSetup, ftime); CHKERRQ(ierr); 1673ccaff030SJeremy L Thompson 1674ccaff030SJeremy L Thompson ierr = VecAXPY(Q, -1.0, Qexact); CHKERRQ(ierr); 1675ccaff030SJeremy L Thompson ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr); 1676cfa64770SLeila Ghaffari CeedVectorDestroy(&q0ceed); 1677ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1678ccaff030SJeremy L Thompson "Max Error: %g\n", 1679ccaff030SJeremy L Thompson (double)norm); CHKERRQ(ierr); 168084d34d69SLeila Ghaffari // Clean up vectors 168184d34d69SLeila Ghaffari ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr); 168284d34d69SLeila Ghaffari ierr = VecDestroy(&Qexact); CHKERRQ(ierr); 1683ccaff030SJeremy L Thompson } 1684ccaff030SJeremy L Thompson 1685ccaff030SJeremy L Thompson // Output Statistics 1686ccaff030SJeremy L Thompson ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr); 1687*dc8efd83SLeila Ghaffari if (!test) { 1688ccaff030SJeremy L Thompson ierr = PetscPrintf(PETSC_COMM_WORLD, 1689ccaff030SJeremy L Thompson "Time integrator took %D time steps to reach final time %g\n", 1690ccaff030SJeremy L Thompson steps, (double)ftime); CHKERRQ(ierr); 1691ccaff030SJeremy L Thompson } 169284d34d69SLeila Ghaffari // Output numerical values from command line 169384d34d69SLeila Ghaffari ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr); 169484d34d69SLeila Ghaffari 1695335ea220SLeila Ghaffari // Compare reference solution values with current test run for CI 1696*dc8efd83SLeila Ghaffari if (test) { 169784d34d69SLeila Ghaffari PetscViewer viewer; 169884d34d69SLeila Ghaffari // Read reference file 169984d34d69SLeila Ghaffari Vec Qref; 170084d34d69SLeila Ghaffari PetscReal error, Qrefnorm; 170184d34d69SLeila Ghaffari ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr); 1702*dc8efd83SLeila Ghaffari ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer); 170384d34d69SLeila Ghaffari CHKERRQ(ierr); 170484d34d69SLeila Ghaffari ierr = VecLoad(Qref, viewer); CHKERRQ(ierr); 170584d34d69SLeila Ghaffari ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 170684d34d69SLeila Ghaffari 170784d34d69SLeila Ghaffari // Compute error with respect to reference solution 170884d34d69SLeila Ghaffari ierr = VecAXPY(Q, -1.0, Qref); CHKERRQ(ierr); 170984d34d69SLeila Ghaffari ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr); 171084d34d69SLeila Ghaffari ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr); 171184d34d69SLeila Ghaffari ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr); 171284d34d69SLeila Ghaffari ierr = VecDestroy(&Qref); CHKERRQ(ierr); 171384d34d69SLeila Ghaffari // Check error 1714*dc8efd83SLeila Ghaffari if (error > testtol) { 171584d34d69SLeila Ghaffari ierr = PetscPrintf(PETSC_COMM_WORLD, 171684d34d69SLeila Ghaffari "Test failed with error norm %g\n", 171784d34d69SLeila Ghaffari (double)error); CHKERRQ(ierr); 171884d34d69SLeila Ghaffari } 171984d34d69SLeila Ghaffari ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); 172084d34d69SLeila Ghaffari } 17219cf88b28Svaleriabarra 1722ccaff030SJeremy L Thompson // Clean up libCEED 17238b982baeSLeila Ghaffari CeedVectorDestroy(&qdata); 1724ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qceed); 1725ccaff030SJeremy L Thompson CeedVectorDestroy(&user->qdotceed); 1726ccaff030SJeremy L Thompson CeedVectorDestroy(&user->gceed); 1727ccaff030SJeremy L Thompson CeedVectorDestroy(&xcorners); 172884d34d69SLeila Ghaffari CeedBasisDestroy(&basisq); 172984d34d69SLeila Ghaffari CeedBasisDestroy(&basisx); 173084d34d69SLeila Ghaffari CeedBasisDestroy(&basisxc); 173184d34d69SLeila Ghaffari CeedElemRestrictionDestroy(&restrictq); 173284d34d69SLeila Ghaffari CeedElemRestrictionDestroy(&restrictx); 173384d34d69SLeila Ghaffari CeedElemRestrictionDestroy(&restrictqdi); 1734ea6e0f84SLeila Ghaffari CeedQFunctionDestroy(&qf_setupVol); 1735ccaff030SJeremy L Thompson CeedQFunctionDestroy(&qf_ics); 1736ea6e0f84SLeila Ghaffari CeedQFunctionDestroy(&qf_rhsVol); 1737ea6e0f84SLeila Ghaffari CeedQFunctionDestroy(&qf_ifunctionVol); 1738777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxSetup); 1739777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxNS); 1740777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxAdvection2d); 1741777ff853SJeremy L Thompson CeedQFunctionContextDestroy(&ctxSurface); 1742ea6e0f84SLeila Ghaffari CeedOperatorDestroy(&op_setupVol); 1743ccaff030SJeremy L Thompson CeedOperatorDestroy(&op_ics); 17446a0edaf9SLeila Ghaffari CeedOperatorDestroy(&user->op_rhs_vol); 17456a0edaf9SLeila Ghaffari CeedOperatorDestroy(&user->op_ifunction_vol); 17466a0edaf9SLeila Ghaffari CeedDestroy(&ceed); 17476a0edaf9SLeila Ghaffari CeedBasisDestroy(&basisqSur); 17486a0edaf9SLeila Ghaffari CeedBasisDestroy(&basisxSur); 17496a0edaf9SLeila Ghaffari CeedBasisDestroy(&basisxcSur); 17506a0edaf9SLeila Ghaffari CeedQFunctionDestroy(&qf_setupSur); 17517ed8b9c7SLeila Ghaffari CeedQFunctionDestroy(&qf_applySur); 1752ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_rhs); 1753ccaff030SJeremy L Thompson CeedOperatorDestroy(&user->op_ifunction); 1754ccaff030SJeremy L Thompson 1755ccaff030SJeremy L Thompson // Clean up PETSc 1756ccaff030SJeremy L Thompson ierr = VecDestroy(&Q); CHKERRQ(ierr); 1757ccaff030SJeremy L Thompson ierr = VecDestroy(&user->M); CHKERRQ(ierr); 1758ccaff030SJeremy L Thompson ierr = MatDestroy(&interpviz); CHKERRQ(ierr); 1759ccaff030SJeremy L Thompson ierr = DMDestroy(&dmviz); CHKERRQ(ierr); 1760ccaff030SJeremy L Thompson ierr = TSDestroy(&ts); CHKERRQ(ierr); 1761ccaff030SJeremy L Thompson ierr = DMDestroy(&dm); CHKERRQ(ierr); 1762ccaff030SJeremy L Thompson ierr = PetscFree(units); CHKERRQ(ierr); 1763ccaff030SJeremy L Thompson ierr = PetscFree(user); CHKERRQ(ierr); 1764ccaff030SJeremy L Thompson return PetscFinalize(); 1765ccaff030SJeremy L Thompson } 1766