xref: /libCEED/examples/fluids/navierstokes.c (revision ebb4b9bddc594471d392c1002224c01788427abd)
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
3284d34d69SLeila Ghaffari //     ./navierstokes -ceed /gpu/occa -problem advection -degree 1
33ccaff030SJeremy L Thompson //
3484d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test explicit -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
3584d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test implicit_stab_none -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
3684d34d69SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test implicit_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. -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 -stab supg
37ccaff030SJeremy L Thompson 
38ccaff030SJeremy L Thompson /// @file
39ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc
40ccaff030SJeremy L Thompson 
41ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n";
42ccaff030SJeremy L Thompson 
43ccaff030SJeremy L Thompson #include <petscts.h>
44ccaff030SJeremy L Thompson #include <petscdmplex.h>
45ccaff030SJeremy L Thompson #include <ceed.h>
46ccaff030SJeremy L Thompson #include <stdbool.h>
47ccaff030SJeremy L Thompson #include <petscsys.h>
48ccaff030SJeremy L Thompson #include "common.h"
49b0137797SLeila Ghaffari #include "setup-boundary.h"
50ccaff030SJeremy L Thompson #include "advection.h"
51ccaff030SJeremy L Thompson #include "advection2d.h"
52ccaff030SJeremy L Thompson #include "densitycurrent.h"
53ccaff030SJeremy L Thompson 
5484d34d69SLeila Ghaffari #if PETSC_VERSION_LT(3,14,0)
5584d34d69SLeila Ghaffari #  define DMPlexGetClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexGetClosureIndices(a,b,c,d,f,g,i)
5684d34d69SLeila Ghaffari #  define DMPlexRestoreClosureIndices(a,b,c,d,e,f,g,h,i) DMPlexRestoreClosureIndices(a,b,c,d,f,g,i)
5784d34d69SLeila Ghaffari #endif
5884d34d69SLeila Ghaffari 
5984d34d69SLeila Ghaffari // MemType Options
6084d34d69SLeila Ghaffari static const char *const memTypes[] = {
6184d34d69SLeila Ghaffari   "host",
6284d34d69SLeila Ghaffari   "device",
6384d34d69SLeila Ghaffari   "memType", "CEED_MEM_", NULL
6484d34d69SLeila Ghaffari };
6584d34d69SLeila Ghaffari 
66ccaff030SJeremy L Thompson // Problem Options
67ccaff030SJeremy L Thompson typedef enum {
68ccaff030SJeremy L Thompson   NS_DENSITY_CURRENT = 0,
69ccaff030SJeremy L Thompson   NS_ADVECTION = 1,
70ccaff030SJeremy L Thompson   NS_ADVECTION2D = 2,
71ccaff030SJeremy L Thompson } problemType;
72ccaff030SJeremy L Thompson static const char *const problemTypes[] = {
73ccaff030SJeremy L Thompson   "density_current",
74ccaff030SJeremy L Thompson   "advection",
75ccaff030SJeremy L Thompson   "advection2d",
7684d34d69SLeila Ghaffari   "problemType", "NS_", NULL
77ccaff030SJeremy L Thompson };
78ccaff030SJeremy L Thompson 
791184866aSLeila Ghaffari // Wind Options for Advection
801184866aSLeila Ghaffari typedef enum {
811184866aSLeila Ghaffari   ADVECTION_WIND_ROTATION = 0,
821184866aSLeila Ghaffari   ADVECTION_WIND_TRANSLATION = 1,
831184866aSLeila Ghaffari } WindType;
841184866aSLeila Ghaffari static const char *const WindTypes[] = {
851184866aSLeila Ghaffari   "rotation",
861184866aSLeila Ghaffari   "translation",
871184866aSLeila Ghaffari   "WindType", "ADVECTION_WIND_", NULL
881184866aSLeila Ghaffari };
891184866aSLeila Ghaffari 
90ccaff030SJeremy L Thompson typedef enum {
91ccaff030SJeremy L Thompson   STAB_NONE = 0,
92ccaff030SJeremy L Thompson   STAB_SU = 1,   // Streamline Upwind
93ccaff030SJeremy L Thompson   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
94ccaff030SJeremy L Thompson } StabilizationType;
95ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = {
9684d34d69SLeila Ghaffari   "none",
97ccaff030SJeremy L Thompson   "SU",
98ccaff030SJeremy L Thompson   "SUPG",
99ccaff030SJeremy L Thompson   "StabilizationType", "STAB_", NULL
100ccaff030SJeremy L Thompson };
101ccaff030SJeremy L Thompson 
10284d34d69SLeila Ghaffari // Test Options
10384d34d69SLeila Ghaffari typedef enum {
10484d34d69SLeila Ghaffari   TEST_NONE = 0,               // Non test mode
10584d34d69SLeila Ghaffari   TEST_EXPLICIT = 1,           // Explicit test
10684d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_NONE = 2, // Implicit test no stab
10784d34d69SLeila Ghaffari   TEST_IMPLICIT_STAB_SUPG = 3, // Implicit test supg stab
10884d34d69SLeila Ghaffari } testType;
10984d34d69SLeila Ghaffari static const char *const testTypes[] = {
11084d34d69SLeila Ghaffari   "none",
11184d34d69SLeila Ghaffari   "explicit",
11284d34d69SLeila Ghaffari   "implicit_stab_none",
11384d34d69SLeila Ghaffari   "implicit_stab_supg",
11484d34d69SLeila Ghaffari   "testType", "TEST_", NULL
11584d34d69SLeila Ghaffari };
11684d34d69SLeila Ghaffari 
11784d34d69SLeila Ghaffari // Tests specific data
11884d34d69SLeila Ghaffari typedef struct {
11984d34d69SLeila Ghaffari   PetscScalar testtol;
12084d34d69SLeila Ghaffari   const char *filepath;
12184d34d69SLeila Ghaffari } testData;
12284d34d69SLeila Ghaffari 
12384d34d69SLeila Ghaffari testData testOptions[] = {
12484d34d69SLeila Ghaffari   [TEST_NONE] = {
12584d34d69SLeila Ghaffari     .testtol = 0.,
12684d34d69SLeila Ghaffari     .filepath = NULL
12784d34d69SLeila Ghaffari   },
12884d34d69SLeila Ghaffari   [TEST_EXPLICIT] = {
12984d34d69SLeila Ghaffari     .testtol = 1E-5,
13084d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-explicit.bin"
13184d34d69SLeila Ghaffari   },
13284d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_NONE] = {
13384d34d69SLeila Ghaffari     .testtol = 5E-4,
13484d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-none.bin"
13584d34d69SLeila Ghaffari   },
13684d34d69SLeila Ghaffari   [TEST_IMPLICIT_STAB_SUPG] = {
13784d34d69SLeila Ghaffari     .testtol = 5E-4,
13884d34d69SLeila Ghaffari     .filepath = "examples/fluids/tests-output/fluids-navierstokes-implicit-stab-supg.bin"
13984d34d69SLeila Ghaffari   }
14084d34d69SLeila Ghaffari };
14184d34d69SLeila Ghaffari 
142ccaff030SJeremy L Thompson // Problem specific data
143ccaff030SJeremy L Thompson typedef struct {
1448b982baeSLeila Ghaffari   CeedInt dim, qdatasizeVol, qdatasizeSur;
145*ebb4b9bdSLeila Ghaffari   CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applyVol_ifunction,
146*ebb4b9bdSLeila Ghaffari                     applySur;
147ccaff030SJeremy L Thompson   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
148ccaff030SJeremy L Thompson                        PetscScalar[], void *);
1497659d40cSLeila Ghaffari   const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc,
1507659d40cSLeila Ghaffari         *applyVol_ifunction_loc, *applySur_loc;
151ccaff030SJeremy L Thompson   const bool non_zero_time;
152ccaff030SJeremy L Thompson } problemData;
153ccaff030SJeremy L Thompson 
154ccaff030SJeremy L Thompson problemData problemOptions[] = {
155ccaff030SJeremy L Thompson   [NS_DENSITY_CURRENT] = {
156ccaff030SJeremy L Thompson     .dim                       = 3,
157ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1588b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
159b0137797SLeila Ghaffari     .setupVol                  = Setup,
160b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
161356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
162356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
163ccaff030SJeremy L Thompson     .ics                       = ICsDC,
164ccaff030SJeremy L Thompson     .ics_loc                   = ICsDC_loc,
165c96c872fSLeila Ghaffari     .applyVol_rhs              = DC,
166c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = DC_loc,
167c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_DC,
168c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_DC_loc,
169ccaff030SJeremy L Thompson     .bc                        = Exact_DC,
17084d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
171ccaff030SJeremy L Thompson   },
172ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
173ccaff030SJeremy L Thompson     .dim                       = 3,
174ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
1758b982baeSLeila Ghaffari     .qdatasizeSur              = 4,
176b0137797SLeila Ghaffari     .setupVol                  = Setup,
177b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
178356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
179356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
180ccaff030SJeremy L Thompson     .ics                       = ICsAdvection,
181ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection_loc,
182c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection,
183c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection_loc,
184c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection,
185c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection_loc,
1867659d40cSLeila Ghaffari     .applySur                  = Advection_Sur,
1877659d40cSLeila Ghaffari     .applySur_loc              = Advection_Sur_loc,
188ccaff030SJeremy L Thompson     .bc                        = Exact_Advection,
18984d34d69SLeila Ghaffari     .non_zero_time             = PETSC_FALSE,
190ccaff030SJeremy L Thompson   },
191ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
192ccaff030SJeremy L Thompson     .dim                       = 2,
193ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 5,
1948b982baeSLeila Ghaffari     .qdatasizeSur              = 3,
195c96c872fSLeila Ghaffari     .setupVol                  = Setup2d,
196c96c872fSLeila Ghaffari     .setupVol_loc              = Setup2d_loc,
197b0137797SLeila Ghaffari     .setupSur                  = SetupBoundary2d,
198b0137797SLeila Ghaffari     .setupSur_loc              = SetupBoundary2d_loc,
199ccaff030SJeremy L Thompson     .ics                       = ICsAdvection2d,
200ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection2d_loc,
201c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection2d,
202c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection2d_loc,
203c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection2d,
204c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection2d_loc,
2057659d40cSLeila Ghaffari     .applySur                  = Advection2d_Sur,
2067659d40cSLeila Ghaffari     .applySur_loc              = Advection2d_Sur_loc,
207ccaff030SJeremy L Thompson     .bc                        = Exact_Advection2d,
20884d34d69SLeila Ghaffari     .non_zero_time             = PETSC_TRUE,
209ccaff030SJeremy L Thompson   },
210ccaff030SJeremy L Thompson };
211ccaff030SJeremy L Thompson 
212ccaff030SJeremy L Thompson // PETSc user data
213ccaff030SJeremy L Thompson typedef struct User_ *User;
214ccaff030SJeremy L Thompson typedef struct Units_ *Units;
215ccaff030SJeremy L Thompson 
216ccaff030SJeremy L Thompson struct User_ {
217ccaff030SJeremy L Thompson   MPI_Comm comm;
218ccaff030SJeremy L Thompson   PetscInt outputfreq;
219ccaff030SJeremy L Thompson   DM dm;
220ccaff030SJeremy L Thompson   DM dmviz;
221ccaff030SJeremy L Thompson   Mat interpviz;
222ccaff030SJeremy L Thompson   Ceed ceed;
223ccaff030SJeremy L Thompson   Units units;
224ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
2251e150236SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction;
226ccaff030SJeremy L Thompson   Vec M;
227ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
228ccaff030SJeremy L Thompson   PetscInt contsteps;
229ccaff030SJeremy L Thompson };
230ccaff030SJeremy L Thompson 
231ccaff030SJeremy L Thompson struct Units_ {
232ccaff030SJeremy L Thompson   // fundamental units
233ccaff030SJeremy L Thompson   PetscScalar meter;
234ccaff030SJeremy L Thompson   PetscScalar kilogram;
235ccaff030SJeremy L Thompson   PetscScalar second;
236ccaff030SJeremy L Thompson   PetscScalar Kelvin;
237ccaff030SJeremy L Thompson   // derived units
238ccaff030SJeremy L Thompson   PetscScalar Pascal;
239ccaff030SJeremy L Thompson   PetscScalar JperkgK;
240ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
241ccaff030SJeremy L Thompson   PetscScalar WpermK;
242ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
243ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
244ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
24516c0476cSLeila Ghaffari   PetscScalar Joule;
246ccaff030SJeremy L Thompson };
247ccaff030SJeremy L Thompson 
248ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
249ccaff030SJeremy L Thompson struct SimpleBC_ {
2507659d40cSLeila Ghaffari   PetscInt nwall, nslip[3];
2517659d40cSLeila Ghaffari   PetscInt walls[6], slips[3][6];
25284d34d69SLeila Ghaffari   PetscBool userbc;
253ccaff030SJeremy L Thompson };
254ccaff030SJeremy L Thompson 
255ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
256ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
257ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
258ccaff030SJeremy L Thompson }
259ccaff030SJeremy L Thompson 
260ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
261ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
26284d34d69SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value,
263ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
264ccaff030SJeremy L Thompson 
265ccaff030SJeremy L Thompson   PetscSection section;
2661184866aSLeila Ghaffari   PetscInt p, Nelem, Ndof, *erestrict, eoffset, nfields, dim, depth;
2670c6c0b13SLeila Ghaffari   DMLabel depthLabel;
2680c6c0b13SLeila Ghaffari   IS depthIS, iterIS;
26984d34d69SLeila Ghaffari   Vec Uloc;
2700c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
271ccaff030SJeremy L Thompson   PetscErrorCode ierr;
272ccaff030SJeremy L Thompson 
273ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
274ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
275da51bdd9SLeila Ghaffari   dim -= height;
276ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
277ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
278ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
279ccaff030SJeremy L Thompson   fieldoff[0] = 0;
280ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
281ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
282ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
283ccaff030SJeremy L Thompson   }
284ccaff030SJeremy L Thompson 
2850c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2860c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2870c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2880c6c0b13SLeila Ghaffari   if (domainLabel) {
2890c6c0b13SLeila Ghaffari     IS domainIS;
2900c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2910c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2920c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2930c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2940c6c0b13SLeila Ghaffari   } else {
2950c6c0b13SLeila Ghaffari     iterIS = depthIS;
2960c6c0b13SLeila Ghaffari   }
2970c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2980c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
299ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
3000c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
3010c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
302ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
30384d34d69SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, PETSC_TRUE,
30484d34d69SLeila Ghaffari                                    &numindices, &indices, NULL, NULL);
30584d34d69SLeila Ghaffari     CHKERRQ(ierr);
30632b5ec5fSJed Brown     bool flip = false;
30732b5ec5fSJed Brown     if (height > 0) {
30832b5ec5fSJed Brown       PetscInt numCells, numFaces, start = -1;
30932b5ec5fSJed Brown       const PetscInt *orients, *faces, *cells;
31032b5ec5fSJed Brown       ierr = DMPlexGetSupport(dm, c, &cells); CHKERRQ(ierr);
31132b5ec5fSJed Brown       ierr = DMPlexGetSupportSize(dm, c, &numCells); CHKERRQ(ierr);
312*ebb4b9bdSLeila Ghaffari       if (numCells != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
313*ebb4b9bdSLeila Ghaffari                                     "Expected one cell in support of exterior face, but got %D cells", numCells);
31432b5ec5fSJed Brown       ierr = DMPlexGetCone(dm, cells[0], &faces); CHKERRQ(ierr);
31532b5ec5fSJed Brown       ierr = DMPlexGetConeSize(dm, cells[0], &numFaces); CHKERRQ(ierr);
31632b5ec5fSJed Brown       for (PetscInt i=0; i<numFaces; i++) {if (faces[i] == c) start = i;}
317*ebb4b9bdSLeila Ghaffari       if (start < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT,
318*ebb4b9bdSLeila Ghaffari                                 "Could not find face %D in cone of its support", c);
31932b5ec5fSJed Brown       ierr = DMPlexGetConeOrientation(dm, cells[0], &orients); CHKERRQ(ierr);
32032b5ec5fSJed Brown       if (orients[start] < 0) flip = true;
32132b5ec5fSJed Brown     }
32284d34d69SLeila Ghaffari     if (numindices % fieldoff[nfields]) SETERRQ1(PETSC_COMM_SELF,
32384d34d69SLeila Ghaffari           PETSC_ERR_ARG_INCOMP, "Number of closure indices not compatible with Cell %D",
32484d34d69SLeila Ghaffari           c);
325ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
326ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
32732b5ec5fSJed Brown       PetscInt ii = i;
32832b5ec5fSJed Brown       if (flip) {
32932b5ec5fSJed Brown         if (P == nnodes) ii = nnodes - 1 - i;
33032b5ec5fSJed Brown         else if (P*P == nnodes) {
33132b5ec5fSJed Brown           PetscInt row = i / P, col = i % P;
33232b5ec5fSJed Brown           ii = row + col * P;
333*ebb4b9bdSLeila Ghaffari         } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP,
334*ebb4b9bdSLeila Ghaffari                           "No support for flipping point with %D nodes != P (%D) or P^2", nnodes, P);
33532b5ec5fSJed Brown       }
336ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
337ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
338ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
339ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
34032b5ec5fSJed Brown           if (Involute(indices[fieldoff[f]*nnodes + ii*ncomp[f] + j])
34132b5ec5fSJed Brown               != Involute(indices[ii*ncomp[0]]) + fieldoff[f] + j)
342ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
343ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
34432b5ec5fSJed Brown                      c, ii, f, j);
345ccaff030SJeremy L Thompson         }
346ccaff030SJeremy L Thompson       }
347ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
34832b5ec5fSJed Brown       PetscInt loc = Involute(indices[ii*ncomp[0]]);
3496f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
350ccaff030SJeremy L Thompson     }
35184d34d69SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, PETSC_TRUE,
35284d34d69SLeila Ghaffari                                        &numindices, &indices, NULL, NULL);
35384d34d69SLeila Ghaffari     CHKERRQ(ierr);
354ccaff030SJeremy L Thompson   }
3550c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
3560c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
3570c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
358ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
3590c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
3600c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
3610c6c0b13SLeila Ghaffari 
362ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
363ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
364ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
3656f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
3666f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
3676f55dfd5Svaleriabarra                             Erestrict);
368ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
369ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
370ccaff030SJeremy L Thompson }
371ccaff030SJeremy L Thompson 
372c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
3731e150236SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height,
3741e150236SLeila Ghaffari     DMLabel domainLabel, PetscInt value, CeedInt P, CeedInt Q, CeedInt qdatasize,
3751e150236SLeila Ghaffari     CeedElemRestriction *restrictq, CeedElemRestriction *restrictx,
3761e150236SLeila Ghaffari     CeedElemRestriction *restrictqdi) {
377c96c872fSLeila Ghaffari 
378c96c872fSLeila Ghaffari   DM dmcoord;
3791e150236SLeila Ghaffari   CeedInt dim, localNelem;
3801e150236SLeila Ghaffari   CeedInt Qdim;
381c96c872fSLeila Ghaffari   PetscErrorCode ierr;
382c96c872fSLeila Ghaffari 
383c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
3841e150236SLeila Ghaffari   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
3851e150236SLeila Ghaffari   dim -= height;
3861e150236SLeila Ghaffari   Qdim = CeedIntPow(Q, dim);
387c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
388c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
389c96c872fSLeila Ghaffari   CHKERRQ(ierr);
390*ebb4b9bdSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value,
391*ebb4b9bdSLeila Ghaffari                                    restrictq);
392c96c872fSLeila Ghaffari   CHKERRQ(ierr);
393*ebb4b9bdSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value,
394*ebb4b9bdSLeila Ghaffari                                    restrictx);
395c96c872fSLeila Ghaffari   CHKERRQ(ierr);
396c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
397c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
398c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
399c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
400c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
401c96c872fSLeila Ghaffari }
402c96c872fSLeila Ghaffari 
4031e150236SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
4047659d40cSLeila Ghaffari static PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc,
4057659d40cSLeila Ghaffari     WindType wind_type, CeedOperator op_applyVol, CeedQFunction qf_applySur,
4067659d40cSLeila Ghaffari     CeedQFunction qf_setupSur, CeedInt height, CeedInt numP_Sur, CeedInt numQ_Sur,
4077659d40cSLeila Ghaffari     CeedInt qdatasizeSur, CeedInt NqptsSur, CeedBasis basisxSur,
4087659d40cSLeila Ghaffari     CeedBasis basisqSur, CeedOperator *op_apply) {
409ca3ac6ddSLeila Ghaffari 
4107659d40cSLeila Ghaffari   CeedInt dim, nFace;
4117659d40cSLeila Ghaffari   PetscInt lsize, localNelemSur[6];
4121e150236SLeila Ghaffari   Vec Xloc;
4137659d40cSLeila Ghaffari   CeedVector xcorners, qdataSur[6];
4147659d40cSLeila Ghaffari   CeedOperator op_setupSur[6], op_applySur[6];
4157659d40cSLeila Ghaffari   CeedElemRestriction restrictxSur[6], restrictqSur[6], restrictqdiSur[6];
416ca3ac6ddSLeila Ghaffari   DMLabel domainLabel;
4171e150236SLeila Ghaffari   PetscScalar *x;
418ca3ac6ddSLeila Ghaffari   PetscErrorCode ierr;
419ca3ac6ddSLeila Ghaffari 
420ca3ac6ddSLeila Ghaffari   PetscFunctionBeginUser;
421ca3ac6ddSLeila Ghaffari   // Composite Operaters
422ca3ac6ddSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
423*ebb4b9bdSLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply,
424*ebb4b9bdSLeila Ghaffari                               op_applyVol); // Apply a Sub-Operator for the volume
425ca3ac6ddSLeila Ghaffari 
4267659d40cSLeila Ghaffari   if (wind_type == ADVECTION_WIND_TRANSLATION) {
4277659d40cSLeila Ghaffari     bc->nwall = 0;
4287659d40cSLeila Ghaffari     bc->nslip[0] = bc->nslip[1] = bc->nslip[2] = 0;
4291e150236SLeila Ghaffari     ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
4301e150236SLeila Ghaffari     ierr = VecGetLocalSize(Xloc, &lsize); CHKERRQ(ierr);
4311e150236SLeila Ghaffari     ierr = CeedVectorCreate(ceed, lsize, &xcorners); CHKERRQ(ierr);
4321e150236SLeila Ghaffari     ierr = VecGetArray(Xloc, &x); CHKERRQ(ierr);
4331e150236SLeila Ghaffari     CeedVectorSetArray(xcorners, CEED_MEM_HOST, CEED_USE_POINTER, x);
434ca3ac6ddSLeila Ghaffari     ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
4357659d40cSLeila Ghaffari     ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
4367659d40cSLeila Ghaffari     if (dim == 2) nFace = 4;
4377659d40cSLeila Ghaffari     if (dim == 3) nFace = 6;
4389fe13df9SLeila Ghaffari 
4397659d40cSLeila Ghaffari     // Create CEED Operator for each boundary face
4407659d40cSLeila Ghaffari     for (CeedInt i=0; i<nFace; i++) {
4417659d40cSLeila Ghaffari       ierr = GetRestrictionForDomain(ceed, dm, height, domainLabel, i+1, numP_Sur,
4427659d40cSLeila Ghaffari                                      numQ_Sur, qdatasizeSur, &restrictqSur[i], &restrictxSur[i],
4437659d40cSLeila Ghaffari                                      &restrictqdiSur[i]); CHKERRQ(ierr);
4447659d40cSLeila Ghaffari       // Create the CEED vectors that will be needed in Boundary setup
4457659d40cSLeila Ghaffari       CeedElemRestrictionGetNumElements(restrictqSur[i], &localNelemSur[i]);
4467659d40cSLeila Ghaffari       CeedVectorCreate(ceed, qdatasizeSur*localNelemSur[i]*NqptsSur, &qdataSur[i]);
4477659d40cSLeila Ghaffari       // Create the operator that builds the quadrature data for the Boundary operator
4487659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur[i]);
449*ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "dx", restrictxSur[i], basisxSur,
450*ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4517659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "weight", CEED_ELEMRESTRICTION_NONE,
452ca3ac6ddSLeila Ghaffari                            basisxSur, CEED_VECTOR_NONE);
4537659d40cSLeila Ghaffari       CeedOperatorSetField(op_setupSur[i], "qdataSur", restrictqdiSur[i],
454ca3ac6ddSLeila Ghaffari                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
4557659d40cSLeila Ghaffari       // Create Boundary operator
4567659d40cSLeila Ghaffari       CeedOperatorCreate(ceed, qf_applySur, NULL, NULL, &op_applySur[i]);
457*ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "q", restrictqSur[i], basisqSur,
458*ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4597659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "qdataSur", restrictqdiSur[i],
4607659d40cSLeila Ghaffari                            CEED_BASIS_COLLOCATED, qdataSur[i]);
4617659d40cSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "x", restrictxSur[i], basisxSur, xcorners);
462*ebb4b9bdSLeila Ghaffari       CeedOperatorSetField(op_applySur[i], "v", restrictqSur[i], basisqSur,
463*ebb4b9bdSLeila Ghaffari                            CEED_VECTOR_ACTIVE);
4647659d40cSLeila Ghaffari       // Apply CEED operator for Boundary setup
465*ebb4b9bdSLeila Ghaffari       CeedOperatorApply(op_setupSur[i], xcorners, qdataSur[i],
466*ebb4b9bdSLeila Ghaffari                         CEED_REQUEST_IMMEDIATE);
4677659d40cSLeila Ghaffari       // Apply Sub-Operator for the Boundary
4687659d40cSLeila Ghaffari       CeedCompositeOperatorAddSub(*op_apply, op_applySur[i]);
4699fe13df9SLeila Ghaffari     }
4701e150236SLeila Ghaffari     CeedVectorDestroy(&xcorners);
471ca3ac6ddSLeila Ghaffari   }
472ca3ac6ddSLeila Ghaffari   PetscFunctionReturn(0);
473ca3ac6ddSLeila Ghaffari }
474ca3ac6ddSLeila Ghaffari 
475ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
476ccaff030SJeremy L Thompson   PetscErrorCode ierr;
477ccaff030SJeremy L Thompson   PetscInt m;
478ccaff030SJeremy L Thompson 
479ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
480ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
481ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
482ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
483ccaff030SJeremy L Thompson }
484ccaff030SJeremy L Thompson 
485ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
486ccaff030SJeremy L Thompson   PetscErrorCode ierr;
487ccaff030SJeremy L Thompson   PetscInt mceed, mpetsc;
488ccaff030SJeremy L Thompson   PetscScalar *a;
489ccaff030SJeremy L Thompson 
490ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
491ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
492ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
493ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
49484d34d69SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",
49584d34d69SLeila Ghaffari                                   mpetsc, mceed);
496ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
497ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
498ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
499ccaff030SJeremy L Thompson }
500ccaff030SJeremy L Thompson 
501ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
502ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
503ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
504ccaff030SJeremy L Thompson   PetscErrorCode ierr;
505ccaff030SJeremy L Thompson   Vec Qbc;
506ccaff030SJeremy L Thompson 
507ccaff030SJeremy L Thompson   PetscFunctionBegin;
508ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
509ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
510ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
511ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
512ccaff030SJeremy L Thompson }
513ccaff030SJeremy L Thompson 
514ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
515ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
516ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
517ccaff030SJeremy L Thompson   PetscErrorCode ierr;
518ccaff030SJeremy L Thompson   User user = *(User *)userData;
519ccaff030SJeremy L Thompson   PetscScalar *q, *g;
520ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
521ccaff030SJeremy L Thompson 
522ccaff030SJeremy L Thompson   // Global-to-local
523ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
524ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
525ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
526ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
527ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
528ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
529ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
530ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
531ccaff030SJeremy L Thompson 
532ccaff030SJeremy L Thompson   // Ceed Vectors
533ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
534ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
535ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
536ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
537ccaff030SJeremy L Thompson 
538ccaff030SJeremy L Thompson   // Apply CEED operator
539ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
540ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
541ccaff030SJeremy L Thompson 
542ccaff030SJeremy L Thompson   // Restore vectors
543ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
544ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
545ccaff030SJeremy L Thompson 
546ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
547ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
548ccaff030SJeremy L Thompson 
549ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
550ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
551ccaff030SJeremy L Thompson   CHKERRQ(ierr);
552ccaff030SJeremy L Thompson 
553ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
554ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
555ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
556ccaff030SJeremy L Thompson }
557ccaff030SJeremy L Thompson 
558ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
559ccaff030SJeremy L Thompson                                    void *userData) {
560ccaff030SJeremy L Thompson   PetscErrorCode ierr;
561ccaff030SJeremy L Thompson   User user = *(User *)userData;
562ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
563ccaff030SJeremy L Thompson   PetscScalar *g;
564ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
565ccaff030SJeremy L Thompson 
566ccaff030SJeremy L Thompson   // Global-to-local
567ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
568ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
569ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
570ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
571ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
572ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
573ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
574ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
575ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
576ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
577ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
578ccaff030SJeremy L Thompson 
579ccaff030SJeremy L Thompson   // Ceed Vectors
580ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
581ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
582ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
583ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
584ccaff030SJeremy L Thompson                      (PetscScalar *)q);
585ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
586ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
587ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
588ccaff030SJeremy L Thompson 
589ccaff030SJeremy L Thompson   // Apply CEED operator
590ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
591ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
592ccaff030SJeremy L Thompson 
593ccaff030SJeremy L Thompson   // Restore vectors
594ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
595ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
596ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
597ccaff030SJeremy L Thompson 
598ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
599ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
600ccaff030SJeremy L Thompson 
601ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
602ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
603ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
604ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
605ccaff030SJeremy L Thompson }
606ccaff030SJeremy L Thompson 
607ccaff030SJeremy L Thompson // User provided TS Monitor
608ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
609ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
610ccaff030SJeremy L Thompson   User user = ctx;
611ccaff030SJeremy L Thompson   Vec Qloc;
612ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
613ccaff030SJeremy L Thompson   PetscViewer viewer;
614ccaff030SJeremy L Thompson   PetscErrorCode ierr;
615ccaff030SJeremy L Thompson 
616ccaff030SJeremy L Thompson   // Set up output
617ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
618ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
619ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
620ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
621ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
622ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
623ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
624ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
625ccaff030SJeremy L Thompson 
626ccaff030SJeremy L Thompson   // Output
627ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
628ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
629ccaff030SJeremy L Thompson   CHKERRQ(ierr);
630ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
631ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
632ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
6339d801c56SJed Brown   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
634ccaff030SJeremy L Thompson   if (user->dmviz) {
635ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
636ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
637ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
638ccaff030SJeremy L Thompson 
639ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
640ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
641ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
642ccaff030SJeremy L Thompson     CHKERRQ(ierr);
643ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
644ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
645ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
646ccaff030SJeremy L Thompson     CHKERRQ(ierr);
647ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
648ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
649ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
650ccaff030SJeremy L Thompson     CHKERRQ(ierr);
651ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
652ccaff030SJeremy L Thompson                               filepath_refined,
653ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
654ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
655ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
656ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
657ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
658ccaff030SJeremy L Thompson   }
659ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
660ccaff030SJeremy L Thompson 
661ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
662ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
663ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
664ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
665ccaff030SJeremy L Thompson   CHKERRQ(ierr);
666ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
667ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
668ccaff030SJeremy L Thompson 
669ccaff030SJeremy L Thompson   // Save time stamp
670ccaff030SJeremy L Thompson   // Dimensionalize time back
671ccaff030SJeremy L Thompson   time /= user->units->second;
672ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
673ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
674ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
675ccaff030SJeremy L Thompson   CHKERRQ(ierr);
676ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
677ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
678ccaff030SJeremy L Thompson   #else
679ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
680ccaff030SJeremy L Thompson   #endif
681ccaff030SJeremy L Thompson   CHKERRQ(ierr);
682ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
683ccaff030SJeremy L Thompson 
684ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
685ccaff030SJeremy L Thompson }
686ccaff030SJeremy L Thompson 
68784d34d69SLeila Ghaffari static PetscErrorCode ICs_FixMultiplicity(CeedOperator op_ics,
688ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
689ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
690ccaff030SJeremy L Thompson   PetscErrorCode ierr;
691ccaff030SJeremy L Thompson   CeedVector multlvec;
692ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
693ccaff030SJeremy L Thompson 
694ccaff030SJeremy L Thompson   ctxSetup->time = time;
695ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
696ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
697ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
698ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
699ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
700ccaff030SJeremy L Thompson 
701ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
702ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
703ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
704ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
705ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
706ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
707ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
708ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
709ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
710ccaff030SJeremy L Thompson   CHKERRQ(ierr);
711ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
712ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
713ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
714ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
715ccaff030SJeremy L Thompson 
716ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
717ccaff030SJeremy L Thompson }
718ccaff030SJeremy L Thompson 
719ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
720ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
721ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
722ccaff030SJeremy L Thompson   PetscErrorCode ierr;
723ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
724ccaff030SJeremy L Thompson   CeedOperator op_mass;
725ccaff030SJeremy L Thompson   CeedVector mceed;
726ccaff030SJeremy L Thompson   Vec Mloc;
727ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
728ccaff030SJeremy L Thompson 
729ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
730ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
731ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
732ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
733ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
734ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
735ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
736ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
737ccaff030SJeremy L Thompson 
738ccaff030SJeremy L Thompson   // Create the mass operator
739ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
740ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
741ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
742ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
743ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
744ccaff030SJeremy L Thompson 
745ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
746ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
747ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
748ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
749ccaff030SJeremy L Thompson 
750ccaff030SJeremy L Thompson   {
751ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
752ccaff030SJeremy L Thompson     CeedVector onesvec;
753ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
754ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
755ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
756ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
757ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
758ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
759ccaff030SJeremy L Thompson   }
760ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
761ccaff030SJeremy L Thompson 
762ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
763ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
764ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
765ccaff030SJeremy L Thompson 
766ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
767ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
768ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
769ccaff030SJeremy L Thompson }
770ccaff030SJeremy L Thompson 
77184d34d69SLeila Ghaffari static PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
772ff6701fcSJed Brown                               SimpleBC bc, void *ctxSetup) {
773ccaff030SJeremy L Thompson   PetscErrorCode ierr;
774ccaff030SJeremy L Thompson 
775ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
776ccaff030SJeremy L Thompson   {
777ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
778ccaff030SJeremy L Thompson     PetscFE fe;
779ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
780ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
781ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
78232ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
783ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
784ccaff030SJeremy L Thompson     ierr = DMAddField(dm, NULL,(PetscObject)fe); CHKERRQ(ierr);
785ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
78607af6069Svaleriabarra     {
78707af6069Svaleriabarra       PetscInt comps[1] = {1};
78807af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
78907af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
79007af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
79107af6069Svaleriabarra       comps[0] = 2;
79207af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
79307af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
79407af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
79507af6069Svaleriabarra       comps[0] = 3;
79607af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
79707af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
79807af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
79907af6069Svaleriabarra     }
80084d34d69SLeila Ghaffari     if (bc->userbc == PETSC_TRUE) {
80184d34d69SLeila Ghaffari       for (PetscInt c = 0; c < 3; c++) {
80284d34d69SLeila Ghaffari         for (PetscInt s = 0; s < bc->nslip[c]; s++) {
80384d34d69SLeila Ghaffari           for (PetscInt w = 0; w < bc->nwall; w++) {
80484d34d69SLeila Ghaffari             if (bc->slips[c][s] == bc->walls[w])
80584d34d69SLeila Ghaffari               SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
80684d34d69SLeila Ghaffari                        "Boundary condition already set on face %D!\n", bc->walls[w]);
80784d34d69SLeila Ghaffari 
80884d34d69SLeila Ghaffari           }
80984d34d69SLeila Ghaffari         }
81084d34d69SLeila Ghaffari       }
81184d34d69SLeila Ghaffari     }
81284d34d69SLeila Ghaffari     // Wall boundary conditions are zero energy density and zero flux for
81384d34d69SLeila Ghaffari     //   velocity in advection/advection2d, and zero velocity and zero flux
81484d34d69SLeila Ghaffari     //   for mass density and energy density in density_current
81584d34d69SLeila Ghaffari     {
81684d34d69SLeila Ghaffari       if (problem->bc == Exact_Advection || problem->bc == Exact_Advection2d) {
81784d34d69SLeila Ghaffari         PetscInt comps[1] = {4};
81884d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
81984d34d69SLeila Ghaffari                              1, comps, (void(*)(void))problem->bc,
82084d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
82184d34d69SLeila Ghaffari       } else if (problem->bc == Exact_DC) {
82284d34d69SLeila Ghaffari         PetscInt comps[3] = {1, 2, 3};
82384d34d69SLeila Ghaffari         ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
82484d34d69SLeila Ghaffari                              3, comps, (void(*)(void))problem->bc,
82584d34d69SLeila Ghaffari                              bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
82684d34d69SLeila Ghaffari       } else
82784d34d69SLeila Ghaffari         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_NULL,
82884d34d69SLeila Ghaffari                 "Undefined boundary conditions for this problem");
82984d34d69SLeila Ghaffari     }
830ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE, NULL);
831ccaff030SJeremy L Thompson     CHKERRQ(ierr);
832ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
833ccaff030SJeremy L Thompson   }
834ccaff030SJeremy L Thompson   {
835ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
836ccaff030SJeremy L Thompson     PetscSection section;
837ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
838ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
839ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
840ccaff030SJeremy L Thompson     CHKERRQ(ierr);
841ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
842ccaff030SJeremy L Thompson     CHKERRQ(ierr);
843ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
844ccaff030SJeremy L Thompson     CHKERRQ(ierr);
845ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
846ccaff030SJeremy L Thompson     CHKERRQ(ierr);
847ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
848ccaff030SJeremy L Thompson     CHKERRQ(ierr);
849ccaff030SJeremy L Thompson   }
850ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
851ccaff030SJeremy L Thompson }
852ccaff030SJeremy L Thompson 
853ccaff030SJeremy L Thompson int main(int argc, char **argv) {
854ccaff030SJeremy L Thompson   PetscInt ierr;
855ccaff030SJeremy L Thompson   MPI_Comm comm;
85684d34d69SLeila Ghaffari   DM dm, dmcoord, dmviz;
857ccaff030SJeremy L Thompson   Mat interpviz;
858ccaff030SJeremy L Thompson   TS ts;
859ccaff030SJeremy L Thompson   TSAdapt adapt;
860ccaff030SJeremy L Thompson   User user;
861ccaff030SJeremy L Thompson   Units units;
862ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
86384d34d69SLeila Ghaffari   PetscInt localNelemVol, lnodes, gnodes, steps;
864ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
865ccaff030SJeremy L Thompson   PetscMPIInt rank;
866ccaff030SJeremy L Thompson   PetscScalar ftime;
867ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
868ccaff030SJeremy L Thompson   Ceed ceed;
86984d34d69SLeila Ghaffari   CeedInt numP, numQ;
870cfa64770SLeila Ghaffari   CeedVector xcorners, qdata, q0ceed;
87184d34d69SLeila Ghaffari   CeedBasis basisx, basisxc, basisq;
87284d34d69SLeila Ghaffari   CeedElemRestriction restrictx, restrictq, restrictqdi;
873cfa64770SLeila Ghaffari   CeedQFunction qf_setupVol, qf_ics, qf_rhsVol, qf_ifunctionVol;
874cfa64770SLeila Ghaffari   CeedOperator op_setupVol, op_ics;
875ccaff030SJeremy L Thompson   CeedScalar Rd;
87684d34d69SLeila Ghaffari   CeedMemType memtyperequested;
877ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
87816c0476cSLeila Ghaffari               kgpersquaredms, Joulepercubicm, Joule;
879ccaff030SJeremy L Thompson   problemType problemChoice;
880ccaff030SJeremy L Thompson   problemData *problem = NULL;
8811184866aSLeila Ghaffari   WindType wind_type;
882ccaff030SJeremy L Thompson   StabilizationType stab;
88384d34d69SLeila Ghaffari   testType testChoice;
88484d34d69SLeila Ghaffari   testData *test = NULL;
88584d34d69SLeila Ghaffari   PetscBool implicit;
886cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
887ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
88884d34d69SLeila Ghaffari     .nslip = {2, 2, 2},
88984d34d69SLeila Ghaffari     .slips = {{5, 6}, {3, 4}, {1, 2}}
890ccaff030SJeremy L Thompson   };
891ccaff030SJeremy L Thompson   double start, cpu_time_used;
89284d34d69SLeila Ghaffari   // Check PETSc CUDA support
89384d34d69SLeila Ghaffari   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
89484d34d69SLeila Ghaffari   // *INDENT-OFF*
89584d34d69SLeila Ghaffari   #ifdef PETSC_HAVE_CUDA
89684d34d69SLeila Ghaffari   petschavecuda = PETSC_TRUE;
89784d34d69SLeila Ghaffari   #else
89884d34d69SLeila Ghaffari   petschavecuda = PETSC_FALSE;
89984d34d69SLeila Ghaffari   #endif
90084d34d69SLeila Ghaffari   // *INDENT-ON*
901ccaff030SJeremy L Thompson 
902ccaff030SJeremy L Thompson   // Create the libCEED contexts
903ccaff030SJeremy L Thompson   PetscScalar meter      = 1e-2;     // 1 meter in scaled length units
904ccaff030SJeremy L Thompson   PetscScalar second     = 1e-2;     // 1 second in scaled time units
905ccaff030SJeremy L Thompson   PetscScalar kilogram   = 1e-6;     // 1 kilogram in scaled mass units
906ccaff030SJeremy L Thompson   PetscScalar Kelvin     = 1;        // 1 Kelvin in scaled temperature units
907ccaff030SJeremy L Thompson   CeedScalar theta0      = 300.;     // K
908ccaff030SJeremy L Thompson   CeedScalar thetaC      = -15.;     // K
909ccaff030SJeremy L Thompson   CeedScalar P0          = 1.e5;     // Pa
91016c0476cSLeila Ghaffari   CeedScalar E_wind      = 1.e6;     // J
911ccaff030SJeremy L Thompson   CeedScalar N           = 0.01;     // 1/s
912ccaff030SJeremy L Thompson   CeedScalar cv          = 717.;     // J/(kg K)
913ccaff030SJeremy L Thompson   CeedScalar cp          = 1004.;    // J/(kg K)
914ccaff030SJeremy L Thompson   CeedScalar g           = 9.81;     // m/s^2
915ccaff030SJeremy L Thompson   CeedScalar lambda      = -2./3.;   // -
916ccaff030SJeremy L Thompson   CeedScalar mu          = 75.;      // Pa s, dynamic viscosity
917ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
918ccaff030SJeremy L Thompson   CeedScalar k           = 0.02638;  // W/(m K)
919ccaff030SJeremy L Thompson   CeedScalar CtauS       = 0.;       // dimensionless
920ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;       // [0,1]
921ccaff030SJeremy L Thompson   PetscScalar lx         = 8000.;    // m
922ccaff030SJeremy L Thompson   PetscScalar ly         = 8000.;    // m
923ccaff030SJeremy L Thompson   PetscScalar lz         = 4000.;    // m
924ccaff030SJeremy L Thompson   CeedScalar rc          = 1000.;    // m (Radius of bubble)
925ccaff030SJeremy L Thompson   PetscScalar resx       = 1000.;    // m (resolution in x)
926ccaff030SJeremy L Thompson   PetscScalar resy       = 1000.;    // m (resolution in y)
927ccaff030SJeremy L Thompson   PetscScalar resz       = 1000.;    // m (resolution in z)
928ccaff030SJeremy L Thompson   PetscInt outputfreq    = 10;       // -
929ccaff030SJeremy L Thompson   PetscInt contsteps     = 0;        // -
93084d34d69SLeila Ghaffari   PetscInt degree        = 1;        // -
93184d34d69SLeila Ghaffari   PetscInt qextra        = 2;        // -
932ea6e0f84SLeila Ghaffari   PetscInt qextraSur     = 2;        // -
9331184866aSLeila Ghaffari   PetscReal center[3], dc_axis[3] = {0, 0, 0}, wind[3] = {1., 0, 0};
934ccaff030SJeremy L Thompson 
935ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
936ccaff030SJeremy L Thompson   if (ierr) return ierr;
937ccaff030SJeremy L Thompson 
938ccaff030SJeremy L Thompson   // Allocate PETSc context
939ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
940ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
941ccaff030SJeremy L Thompson 
942ccaff030SJeremy L Thompson   // Parse command line options
943ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
944ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
945ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
946ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
947ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
948ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
94984d34d69SLeila Ghaffari   testChoice = TEST_NONE;
95084d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-test", "Run tests", NULL,
95184d34d69SLeila Ghaffari                           testTypes, (PetscEnum)testChoice,
95284d34d69SLeila Ghaffari                           (PetscEnum *)&testChoice,
95384d34d69SLeila Ghaffari                           NULL); CHKERRQ(ierr);
95484d34d69SLeila Ghaffari   test = &testOptions[testChoice];
955ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
956ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
957ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
958ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
959ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
9601184866aSLeila Ghaffari   ierr = PetscOptionsEnum("-problem_advection_wind", "Wind type in Advection",
9611184866aSLeila Ghaffari                           NULL, WindTypes, (PetscEnum)(wind_type = ADVECTION_WIND_ROTATION),
9621184866aSLeila Ghaffari                           (PetscEnum *)&wind_type, NULL); CHKERRQ(ierr);
9631184866aSLeila Ghaffari   PetscInt n = problem->dim;
96482c09b01SLeila Ghaffari   PetscBool userWind;
965*ebb4b9bdSLeila Ghaffari   ierr = PetscOptionsRealArray("-problem_advection_wind_translation",
966*ebb4b9bdSLeila Ghaffari                                "Constant wind vector",
96782c09b01SLeila Ghaffari                                NULL, wind, &n, &userWind); CHKERRQ(ierr);
96882c09b01SLeila Ghaffari   if (wind_type == ADVECTION_WIND_ROTATION && userWind) {
969*ebb4b9bdSLeila Ghaffari     ierr = PetscPrintf(comm,
970*ebb4b9bdSLeila Ghaffari                        "Warning! Use -problem_advection_wind_translation only with -problem_advection_wind translation\n");
97182c09b01SLeila Ghaffari     CHKERRQ(ierr);
9721184866aSLeila Ghaffari   }
973*ebb4b9bdSLeila Ghaffari   if (wind_type == ADVECTION_WIND_TRANSLATION
974*ebb4b9bdSLeila Ghaffari       && problemChoice == NS_DENSITY_CURRENT) {
97567babd9cSLeila Ghaffari     SETERRQ(comm, PETSC_ERR_ARG_INCOMP,
97667babd9cSLeila Ghaffari             "-problem_advection_wind translation is not defined for -problem density_current");
97767babd9cSLeila Ghaffari   }
978ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
979ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
980ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
981ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
982ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
983ccaff030SJeremy L Thompson   CHKERRQ(ierr);
98484d34d69SLeila Ghaffari   if (!implicit && stab != STAB_NONE) {
98584d34d69SLeila Ghaffari     ierr = PetscPrintf(comm, "Warning! Use -stab only with -implicit\n");
98684d34d69SLeila Ghaffari     CHKERRQ(ierr);
98784d34d69SLeila Ghaffari   }
988ccaff030SJeremy L Thompson   {
9897573aee6SJed Brown     PetscInt len;
9907573aee6SJed Brown     PetscBool flg;
991ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
992ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
993ccaff030SJeremy L Thompson                                 NULL, bc.walls,
994ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
995ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
9967573aee6SJed Brown     if (flg) {
9977573aee6SJed Brown       bc.nwall = len;
9987573aee6SJed Brown       // Using a no-slip wall disables automatic slip walls (they must be set explicitly)
9997573aee6SJed Brown       bc.nslip[0] = bc.nslip[1] = bc.nslip[2] = 0;
10007573aee6SJed Brown     }
1001ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
1002ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
1003ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
1004ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
1005ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
1006ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
1007ccaff030SJeremy L Thompson                                    &len), &flg);
1008ccaff030SJeremy L Thompson       CHKERRQ(ierr);
100984d34d69SLeila Ghaffari       if (flg) {
101084d34d69SLeila Ghaffari         bc.nslip[j] = len;
101184d34d69SLeila Ghaffari         bc.userbc = PETSC_TRUE;
101284d34d69SLeila Ghaffari       }
1013ccaff030SJeremy L Thompson     }
1014ccaff030SJeremy L Thompson   }
1015cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
1016cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
1017cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
1018ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1019ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
1020ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
1021ccaff030SJeremy L Thompson   meter = fabs(meter);
1022ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
1023ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
1024ccaff030SJeremy L Thompson   second = fabs(second);
1025ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
1026ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
1027ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
1028ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
1029ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
1030ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
1031ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
1032ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
1033ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
1034ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
1035ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
1036ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
1037ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
103816c0476cSLeila Ghaffari   ierr = PetscOptionsScalar("-E_wind", "Total energy of inflow wind",
103916c0476cSLeila Ghaffari                             NULL, E_wind, &E_wind, NULL); CHKERRQ(ierr);
1040ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
1041ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
1042ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
1043ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
1044ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
1045ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
1046ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
1047ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
1048ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
1049ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
1050ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
1051ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
1052ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
1053ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
1054ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
1055ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
1056ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
1057ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
105884d34d69SLeila Ghaffari   if (stab == STAB_NONE && CtauS != 0) {
105984d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
106084d34d69SLeila Ghaffari                        "Warning! Use -CtauS only with -stab su or -stab supg\n");
106184d34d69SLeila Ghaffari     CHKERRQ(ierr);
106284d34d69SLeila Ghaffari   }
1063ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
1064ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
1065ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
1066ccaff030SJeremy L Thompson   CHKERRQ(ierr);
106784d34d69SLeila Ghaffari   if (problemChoice == NS_DENSITY_CURRENT && (CtauS != 0 || strong_form != 0)) {
106884d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
106984d34d69SLeila Ghaffari                        "Warning! Problem density_current does not support -CtauS or -strong_form\n");
107084d34d69SLeila Ghaffari     CHKERRQ(ierr);
107184d34d69SLeila Ghaffari   }
1072ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
1073ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
1074ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
1075ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
1076ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
1077ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
1078ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
1079ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
1080ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
1081ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
1082ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
1083ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
1084ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
1085ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
108682c09b01SLeila Ghaffari   n = problem->dim;
1087ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
1088ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
1089ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
1090ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
1091ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
1092ccaff030SJeremy L Thompson   n = problem->dim;
1093ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
1094ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
1095ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
1096ccaff030SJeremy L Thompson   {
1097ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
1098ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
1099ccaff030SJeremy L Thompson     if (norm > 0) {
1100ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
1101ccaff030SJeremy L Thompson     }
1102ccaff030SJeremy L Thompson   }
1103ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
1104ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
1105ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
1106ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
1107ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
110884d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
110984d34d69SLeila Ghaffari                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
111084d34d69SLeila Ghaffari   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
111184d34d69SLeila Ghaffari                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
111281f92cf0SLeila Ghaffari   PetscBool userQextraSur;
1113*ebb4b9bdSLeila Ghaffari   ierr = PetscOptionsInt("-qextra_boundary",
1114*ebb4b9bdSLeila Ghaffari                          "Number of extra quadrature points on in/outflow faces",
111581f92cf0SLeila Ghaffari                          NULL, qextraSur, &qextraSur, &userQextraSur); CHKERRQ(ierr);
1116*ebb4b9bdSLeila Ghaffari   if ((wind_type == ADVECTION_WIND_ROTATION
1117*ebb4b9bdSLeila Ghaffari        || problemChoice == NS_DENSITY_CURRENT) && userQextraSur) {
1118*ebb4b9bdSLeila Ghaffari     ierr = PetscPrintf(comm,
1119*ebb4b9bdSLeila Ghaffari                        "Warning! Use -qextra_boundary only with -problem_advection_wind translation\n");
112081f92cf0SLeila Ghaffari     CHKERRQ(ierr);
112181f92cf0SLeila Ghaffari   }
112284d34d69SLeila Ghaffari   ierr = PetscStrncpy(user->outputfolder, ".", 2); CHKERRQ(ierr);
1123ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
1124ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
1125ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
112684d34d69SLeila Ghaffari   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
112784d34d69SLeila Ghaffari   ierr = PetscOptionsEnum("-memtype",
112884d34d69SLeila Ghaffari                           "CEED MemType requested", NULL,
112984d34d69SLeila Ghaffari                           memTypes, (PetscEnum)memtyperequested,
113084d34d69SLeila Ghaffari                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
113184d34d69SLeila Ghaffari   CHKERRQ(ierr);
1132ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1133ccaff030SJeremy L Thompson 
1134ccaff030SJeremy L Thompson   // Define derived units
1135ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
1136ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
1137ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
1138ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
1139ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
1140ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
1141ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
114216c0476cSLeila Ghaffari   Joule = kilogram * PetscSqr(meter) / PetscSqr(second);
1143ccaff030SJeremy L Thompson 
1144ccaff030SJeremy L Thompson   // Scale variables to desired units
1145ccaff030SJeremy L Thompson   theta0 *= Kelvin;
1146ccaff030SJeremy L Thompson   thetaC *= Kelvin;
1147ccaff030SJeremy L Thompson   P0 *= Pascal;
114816c0476cSLeila Ghaffari   E_wind *= Joule;
1149ccaff030SJeremy L Thompson   N *= (1./second);
1150ccaff030SJeremy L Thompson   cv *= JperkgK;
1151ccaff030SJeremy L Thompson   cp *= JperkgK;
1152ccaff030SJeremy L Thompson   Rd = cp - cv;
1153ccaff030SJeremy L Thompson   g *= mpersquareds;
1154ccaff030SJeremy L Thompson   mu *= Pascal * second;
1155ccaff030SJeremy L Thompson   k *= WpermK;
1156ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
1157ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
1158ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
1159ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
1160ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
1161ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
1162ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
1163ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
1164ccaff030SJeremy L Thompson 
1165ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
1166cfa64770SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
1167ccaff030SJeremy L Thompson   // Set up the libCEED context
1168ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup = {
1169ccaff030SJeremy L Thompson     .theta0 = theta0,
1170ccaff030SJeremy L Thompson     .thetaC = thetaC,
1171ccaff030SJeremy L Thompson     .P0 = P0,
1172ccaff030SJeremy L Thompson     .N = N,
1173ccaff030SJeremy L Thompson     .cv = cv,
1174ccaff030SJeremy L Thompson     .cp = cp,
1175ccaff030SJeremy L Thompson     .Rd = Rd,
1176ccaff030SJeremy L Thompson     .g = g,
1177ccaff030SJeremy L Thompson     .rc = rc,
1178ccaff030SJeremy L Thompson     .lx = lx,
1179ccaff030SJeremy L Thompson     .ly = ly,
1180ccaff030SJeremy L Thompson     .lz = lz,
1181ccaff030SJeremy L Thompson     .center[0] = center[0],
1182ccaff030SJeremy L Thompson     .center[1] = center[1],
1183ccaff030SJeremy L Thompson     .center[2] = center[2],
1184ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
1185ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
1186ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
11871184866aSLeila Ghaffari     .wind[0] = wind[0],
11881184866aSLeila Ghaffari     .wind[1] = wind[1],
11891184866aSLeila Ghaffari     .wind[2] = wind[2],
1190ccaff030SJeremy L Thompson     .time = 0,
11911184866aSLeila Ghaffari     .wind_type = wind_type,
1192ccaff030SJeremy L Thompson   };
1193ccaff030SJeremy L Thompson 
119484d34d69SLeila Ghaffari   // Create the mesh
1195ccaff030SJeremy L Thompson   {
1196ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
1197ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
119884d34d69SLeila Ghaffari                                NULL, PETSC_TRUE, &dm);
1199ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1200ccaff030SJeremy L Thompson   }
120184d34d69SLeila Ghaffari 
120284d34d69SLeila Ghaffari   // Distribute the mesh over processes
120384d34d69SLeila Ghaffari   {
1204ccaff030SJeremy L Thompson     DM               dmDist = NULL;
1205ccaff030SJeremy L Thompson     PetscPartitioner part;
1206ccaff030SJeremy L Thompson 
1207ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
1208ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1209ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
1210ccaff030SJeremy L Thompson     if (dmDist) {
1211ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
1212ccaff030SJeremy L Thompson       dm  = dmDist;
1213ccaff030SJeremy L Thompson     }
1214ccaff030SJeremy L Thompson   }
1215ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
1216ccaff030SJeremy L Thompson 
121784d34d69SLeila Ghaffari   // Setup DM
1218ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
1219ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
122084d34d69SLeila Ghaffari   ierr = SetUpDM(dm, problem, degree, &bc, &ctxSetup); CHKERRQ(ierr);
122184d34d69SLeila Ghaffari 
122284d34d69SLeila Ghaffari   // Refine DM for high-order viz
1223ccaff030SJeremy L Thompson   dmviz = NULL;
1224ccaff030SJeremy L Thompson   interpviz = NULL;
1225ccaff030SJeremy L Thompson   if (viz_refine) {
1226ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
1227ff6701fcSJed Brown 
1228ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
1229ff6701fcSJed Brown     dmhierarchy[0] = dm;
123084d34d69SLeila Ghaffari     for (PetscInt i = 0, d = degree; i < viz_refine; i++) {
1231ff6701fcSJed Brown       Mat interp_next;
1232ff6701fcSJed Brown 
1233ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
1234ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1235bc6a41f7SJed Brown       ierr = DMClearDS(dmhierarchy[i+1]); CHKERRQ(ierr);
1236bc6a41f7SJed Brown       ierr = DMClearFields(dmhierarchy[i+1]); CHKERRQ(ierr);
1237ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
1238ff6701fcSJed Brown       d = (d + 1) / 2;
1239ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
1240ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
1241ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
1242ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
1243ff6701fcSJed Brown       if (!i) interpviz = interp_next;
1244ff6701fcSJed Brown       else {
1245ff6701fcSJed Brown         Mat C;
1246ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1247ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1248ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1249ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1250ff6701fcSJed Brown         interpviz = C;
1251ff6701fcSJed Brown       }
1252ff6701fcSJed Brown     }
1253cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1254ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1255cb3e2689Svaleriabarra     }
1256ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1257ccaff030SJeremy L Thompson   }
1258ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1259ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1260ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1261ccaff030SJeremy L Thompson   lnodes /= ncompq;
1262ccaff030SJeremy L Thompson 
126384d34d69SLeila Ghaffari   // Initialize CEED
126484d34d69SLeila Ghaffari   CeedInit(ceedresource, &ceed);
126584d34d69SLeila Ghaffari   // Set memtype
126684d34d69SLeila Ghaffari   CeedMemType memtypebackend;
126784d34d69SLeila Ghaffari   CeedGetPreferredMemType(ceed, &memtypebackend);
126884d34d69SLeila Ghaffari   // Check memtype compatibility
126984d34d69SLeila Ghaffari   if (!setmemtyperequest)
127084d34d69SLeila Ghaffari     memtyperequested = memtypebackend;
127184d34d69SLeila Ghaffari   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
127284d34d69SLeila Ghaffari     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
127384d34d69SLeila Ghaffari              "PETSc was not built with CUDA. "
127484d34d69SLeila Ghaffari              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
127584d34d69SLeila Ghaffari 
127684d34d69SLeila Ghaffari   // Set number of 1D nodes and quadrature points
127784d34d69SLeila Ghaffari   numP = degree + 1;
127884d34d69SLeila Ghaffari   numQ = numP + qextra;
127984d34d69SLeila Ghaffari 
128084d34d69SLeila Ghaffari   // Print summary
128184d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1282ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1283ccaff030SJeremy L Thompson     int comm_size;
1284ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1285ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1286ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
128784d34d69SLeila Ghaffari     gnodes = gdofs/ncompq;
1288ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1289ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1290ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
129184d34d69SLeila Ghaffari     const char *usedresource;
129284d34d69SLeila Ghaffari     CeedGetResource(ceed, &usedresource);
1293ccaff030SJeremy L Thompson 
129484d34d69SLeila Ghaffari     ierr = PetscPrintf(comm,
129584d34d69SLeila Ghaffari                        "\n-- Navier-Stokes solver - libCEED + PETSc --\n"
129684d34d69SLeila Ghaffari                        "  rank(s)                              : %d\n"
129784d34d69SLeila Ghaffari                        "  Problem:\n"
129884d34d69SLeila Ghaffari                        "    Problem Name                       : %s\n"
129984d34d69SLeila Ghaffari                        "    Stabilization                      : %s\n"
130084d34d69SLeila Ghaffari                        "  PETSc:\n"
130184d34d69SLeila Ghaffari                        "    Box Faces                          : %s\n"
130284d34d69SLeila Ghaffari                        "  libCEED:\n"
130384d34d69SLeila Ghaffari                        "    libCEED Backend                    : %s\n"
130484d34d69SLeila Ghaffari                        "    libCEED Backend MemType            : %s\n"
130584d34d69SLeila Ghaffari                        "    libCEED User Requested MemType     : %s\n"
130684d34d69SLeila Ghaffari                        "  Mesh:\n"
130784d34d69SLeila Ghaffari                        "    Number of 1D Basis Nodes (P)       : %d\n"
130884d34d69SLeila Ghaffari                        "    Number of 1D Quadrature Points (Q) : %d\n"
130984d34d69SLeila Ghaffari                        "    Global DoFs                        : %D\n"
131084d34d69SLeila Ghaffari                        "    Owned DoFs                         : %D\n"
131184d34d69SLeila Ghaffari                        "    DoFs per node                      : %D\n"
131284d34d69SLeila Ghaffari                        "    Global nodes                       : %D\n"
131384d34d69SLeila Ghaffari                        "    Owned nodes                        : %D\n",
131484d34d69SLeila Ghaffari                        comm_size, problemTypes[problemChoice],
131584d34d69SLeila Ghaffari                        StabilizationTypes[stab], box_faces_str, usedresource,
131684d34d69SLeila Ghaffari                        CeedMemTypes[memtypebackend],
131784d34d69SLeila Ghaffari                        (setmemtyperequest) ?
131884d34d69SLeila Ghaffari                        CeedMemTypes[memtyperequested] : "none",
131984d34d69SLeila Ghaffari                        numP, numQ, gdofs, odofs, ncompq, gnodes, lnodes);
132084d34d69SLeila Ghaffari     CHKERRQ(ierr);
13210c6c0b13SLeila Ghaffari   }
13220c6c0b13SLeila Ghaffari 
1323ccaff030SJeremy L Thompson   // Set up global mass vector
1324ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q, &user->M); CHKERRQ(ierr);
1325ccaff030SJeremy L Thompson 
132684d34d69SLeila Ghaffari   // Set up libCEED
1327ccaff030SJeremy L Thompson   // CEED Bases
1328ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
132984d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP, numQ, CEED_GAUSS,
133084d34d69SLeila Ghaffari                                   &basisq);
133184d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ, CEED_GAUSS,
133284d34d69SLeila Ghaffari                                   &basisx);
133384d34d69SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP,
133484d34d69SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxc);
1335ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1336ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1337ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1338ccaff030SJeremy L Thompson 
1339ccaff030SJeremy L Thompson   // CEED Restrictions
13401e150236SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, numP, numQ,
134184d34d69SLeila Ghaffari                                  qdatasizeVol, &restrictq, &restrictx,
134284d34d69SLeila Ghaffari                                  &restrictqdi); CHKERRQ(ierr);
1343ccaff030SJeremy L Thompson 
1344ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1345ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1346ccaff030SJeremy L Thompson 
1347ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1348bd910870SLeila Ghaffari   CeedInt NqptsVol;
134984d34d69SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisq, &NqptsVol);
135084d34d69SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictq, &localNelemVol);
13518b982baeSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdata);
135284d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &q0ceed, NULL);
1353ccaff030SJeremy L Thompson 
1354ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1355ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1356ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1357ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1358ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
13598b982baeSLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1360ccaff030SJeremy L Thompson 
1361ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1362ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1363ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1364ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1365ccaff030SJeremy L Thompson 
1366ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1367ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1368ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1369ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1370ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1371ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
13728b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1373ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1374ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1375ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1376ccaff030SJeremy L Thompson   }
1377ccaff030SJeremy L Thompson 
1378ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1379ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1380ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1381ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1382ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1383ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1384ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
13858b982baeSLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdata", qdatasizeVol, CEED_EVAL_NONE);
1386ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1387ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1388ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1389ccaff030SJeremy L Thompson   }
1390ccaff030SJeremy L Thompson 
1391ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1392ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
139384d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictx, basisx, CEED_VECTOR_ACTIVE);
1394ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
139584d34d69SLeila Ghaffari                        basisx, CEED_VECTOR_NONE);
139684d34d69SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdata", restrictqdi,
1397ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1398ccaff030SJeremy L Thompson 
1399ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1400ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
140184d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictx, basisxc, CEED_VECTOR_ACTIVE);
140284d34d69SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictq,
1403ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1404ccaff030SJeremy L Thompson 
140584d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qceed, NULL);
140684d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->qdotceed, NULL);
140784d34d69SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictq, &user->gceed, NULL);
1408ccaff030SJeremy L Thompson 
1409ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1410ccaff030SJeremy L Thompson     CeedOperator op;
1411ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
141284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
141384d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
141484d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
14158b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
141684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
141784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
141884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1419d3630711SJed Brown     user->op_rhs_vol = op;
1420ccaff030SJeremy L Thompson   }
1421ccaff030SJeremy L Thompson 
1422ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1423ccaff030SJeremy L Thompson     CeedOperator op;
1424ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
142584d34d69SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
142684d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictq, basisq, CEED_VECTOR_ACTIVE);
142784d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictq, basisq, user->qdotceed);
142884d34d69SLeila Ghaffari     CeedOperatorSetField(op, "qdata", restrictqdi,
14298b982baeSLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdata);
143084d34d69SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictx, basisx, xcorners);
143184d34d69SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
143284d34d69SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictq, basisq, CEED_VECTOR_ACTIVE);
1433d3630711SJed Brown     user->op_ifunction_vol = op;
1434ccaff030SJeremy L Thompson   }
1435ccaff030SJeremy L Thompson 
14366a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
14376a0edaf9SLeila Ghaffari   CeedInt height = 1;
14386a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
14391e150236SLeila Ghaffari   CeedInt numP_Sur = degree + 1;
14401e150236SLeila Ghaffari   CeedInt numQ_Sur = numP_Sur + qextraSur;
1441cfa64770SLeila Ghaffari   const CeedInt qdatasizeSur = problem->qdatasizeSur;
1442cfa64770SLeila Ghaffari   CeedBasis basisxSur, basisxcSur, basisqSur;
1443cfa64770SLeila Ghaffari   CeedInt NqptsSur;
14447ed8b9c7SLeila Ghaffari   CeedQFunction qf_setupSur, qf_applySur;
1445cfa64770SLeila Ghaffari 
1446cfa64770SLeila Ghaffari   // CEED bases for the boundaries
1447*ebb4b9bdSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur,
1448*ebb4b9bdSLeila Ghaffari                                   CEED_GAUSS,
14496a0edaf9SLeila Ghaffari                                   &basisqSur);
14506a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
14516a0edaf9SLeila Ghaffari                                   &basisxSur);
14526a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
14536a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
14546a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
14556a0edaf9SLeila Ghaffari 
1456cfa64770SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the Surface operator
14576a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
14586a0edaf9SLeila Ghaffari                               &qf_setupSur);
14596a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
14606a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
14616a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14626a0edaf9SLeila Ghaffari 
14637659d40cSLeila Ghaffari   // Creat Q-Function for Boundaries
14647ed8b9c7SLeila Ghaffari   qf_applySur = NULL;
14657659d40cSLeila Ghaffari   if (problem->applySur) {
14667659d40cSLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur,
14677ed8b9c7SLeila Ghaffari                                 problem->applySur_loc, &qf_applySur);
14687ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "q", ncompq, CEED_EVAL_INTERP);
14697ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
14707ed8b9c7SLeila Ghaffari     CeedQFunctionAddInput(qf_applySur, "x", ncompx, CEED_EVAL_INTERP);
14717ed8b9c7SLeila Ghaffari     CeedQFunctionAddOutput(qf_applySur, "v", ncompq, CEED_EVAL_INTERP);
14729fe13df9SLeila Ghaffari   }
14739fe13df9SLeila Ghaffari 
14749fe13df9SLeila Ghaffari   // Create CEED Operator for the whole domain
14759fe13df9SLeila Ghaffari   if (!implicit)
1476*ebb4b9bdSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_rhs_vol,
1477*ebb4b9bdSLeila Ghaffari                                    qf_applySur, qf_setupSur,
14787659d40cSLeila Ghaffari                                    height, numP_Sur, numQ_Sur, qdatasizeSur, NqptsSur, basisxSur,
14797659d40cSLeila Ghaffari                                    basisqSur, &user->op_rhs); CHKERRQ(ierr);
14809fe13df9SLeila Ghaffari   if (implicit)
1481*ebb4b9bdSLeila Ghaffari     ierr = CreateOperatorForDomain(ceed, dm, &bc, wind_type, user->op_ifunction_vol,
1482*ebb4b9bdSLeila Ghaffari                                    qf_applySur, qf_setupSur,
14837659d40cSLeila Ghaffari                                    height, numP_Sur, numQ_Sur, qdatasizeSur, NqptsSur, basisxSur,
14847659d40cSLeila Ghaffari                                    basisqSur, &user->op_ifunction); CHKERRQ(ierr);
14857659d40cSLeila Ghaffari   // Set up contex for QFunctions
1486ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1487ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1488ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1489ccaff030SJeremy L Thompson     .CtauS = CtauS,
1490ccaff030SJeremy L Thompson     .strong_form = strong_form,
1491ccaff030SJeremy L Thompson     .stabilization = stab,
1492ccaff030SJeremy L Thompson   };
14937659d40cSLeila Ghaffari   struct SurfaceContext_ ctxSurface = {
149416c0476cSLeila Ghaffari     .E_wind = E_wind,
14957659d40cSLeila Ghaffari     .strong_form = strong_form,
14967659d40cSLeila Ghaffari     .implicit = implicit,
14977659d40cSLeila Ghaffari   };
1498ccaff030SJeremy L Thompson   switch (problemChoice) {
1499ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1500ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1501ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1502ccaff030SJeremy L Thompson           sizeof ctxNS);
1503ccaff030SJeremy L Thompson     break;
1504ccaff030SJeremy L Thompson   case NS_ADVECTION:
1505ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1506ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1507ccaff030SJeremy L Thompson                                              sizeof ctxAdvection2d);
1508ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1509ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1510*ebb4b9bdSLeila Ghaffari     if (qf_applySur) CeedQFunctionSetContext(qf_applySur, &ctxSurface,
1511*ebb4b9bdSLeila Ghaffari           sizeof ctxSurface);
1512ccaff030SJeremy L Thompson   }
1513ccaff030SJeremy L Thompson 
1514ccaff030SJeremy L Thompson   // Set up PETSc context
1515ccaff030SJeremy L Thompson   // Set up units structure
1516ccaff030SJeremy L Thompson   units->meter = meter;
1517ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1518ccaff030SJeremy L Thompson   units->second = second;
1519ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1520ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1521ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1522ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1523ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1524ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1525ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1526ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
152716c0476cSLeila Ghaffari   units->Joule = Joule;
1528ccaff030SJeremy L Thompson 
1529ccaff030SJeremy L Thompson   // Set up user structure
1530ccaff030SJeremy L Thompson   user->comm = comm;
1531ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1532ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1533ccaff030SJeremy L Thompson   user->units = units;
1534ccaff030SJeremy L Thompson   user->dm = dm;
1535ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1536ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1537ccaff030SJeremy L Thompson   user->ceed = ceed;
1538ccaff030SJeremy L Thompson 
15398b982baeSLeila Ghaffari   // Calculate qdata and ICs
1540ccaff030SJeremy L Thompson   // Set up state global and local vectors
1541ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1542ccaff030SJeremy L Thompson 
1543cfa64770SLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
1544ccaff030SJeremy L Thompson 
1545ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1546ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
15478b982baeSLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdata, CEED_REQUEST_IMMEDIATE);
154884d34d69SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictq, basisq, restrictqdi, qdata,
1549ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1550ccaff030SJeremy L Thompson 
155184d34d69SLeila Ghaffari   ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qloc, Q, restrictq,
155284d34d69SLeila Ghaffari                              &ctxSetup, 0.0); CHKERRQ(ierr);
1553ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1554ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1555ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1556ccaff030SJeremy L Thompson     // slower execution.
1557ccaff030SJeremy L Thompson     Vec Qbc;
1558ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1559ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1560ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1561ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1562ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1563ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1564ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
156584d34d69SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_NS);
156684d34d69SLeila Ghaffari     CHKERRQ(ierr);
1567ccaff030SJeremy L Thompson   }
1568ccaff030SJeremy L Thompson 
1569ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1570ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1571ccaff030SJeremy L Thompson   // Gather initial Q values
1572ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1573ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1574ccaff030SJeremy L Thompson     PetscViewer viewer;
1575ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1576ccaff030SJeremy L Thompson     // Read input
1577ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1578ccaff030SJeremy L Thompson                          user->outputfolder);
1579ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1580ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1581ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1582ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1583ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1584ccaff030SJeremy L Thompson   }
1585ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1586ccaff030SJeremy L Thompson 
1587ccaff030SJeremy L Thompson // Create and setup TS
1588ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1589ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1590ccaff030SJeremy L Thompson   if (implicit) {
1591ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1592ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1593ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1594ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1595ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1596ccaff030SJeremy L Thompson     }
1597ccaff030SJeremy L Thompson   } else {
1598ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm, PETSC_ERR_ARG_NULL,
1599ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1600ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1601ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1602ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1603ccaff030SJeremy L Thompson   }
1604ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1605ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1606ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
160784d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {ierr = TSSetMaxSteps(ts, 10); CHKERRQ(ierr);}
1608ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1609ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1610ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1611ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1612ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1613ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
161484d34d69SLeila Ghaffari     if (testChoice == TEST_NONE) {
1615ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1616ccaff030SJeremy L Thompson     }
1617ccaff030SJeremy L Thompson   } else { // continue from time of last output
1618ccaff030SJeremy L Thompson     PetscReal time;
1619ccaff030SJeremy L Thompson     PetscInt count;
1620ccaff030SJeremy L Thompson     PetscViewer viewer;
1621ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1622ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1623ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1624ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1625ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1626ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1627ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1628ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1629ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1630ccaff030SJeremy L Thompson   }
163184d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1632ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1633ccaff030SJeremy L Thompson   }
1634ccaff030SJeremy L Thompson 
1635ccaff030SJeremy L Thompson   // Solve
1636ccaff030SJeremy L Thompson   start = MPI_Wtime();
1637ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1638ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1639ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1640ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts, &ftime); CHKERRQ(ierr);
1641ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1642ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
164384d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1644ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
164584d34d69SLeila Ghaffari                        "Time taken for solution (sec): %g\n",
1646ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1647ccaff030SJeremy L Thompson   }
1648ccaff030SJeremy L Thompson 
1649ccaff030SJeremy L Thompson   // Get error
165084d34d69SLeila Ghaffari   if (problem->non_zero_time && testChoice == TEST_NONE) {
1651ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1652ccaff030SJeremy L Thompson     PetscReal norm;
1653ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1654ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1655ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1656ccaff030SJeremy L Thompson 
165784d34d69SLeila Ghaffari     ierr = ICs_FixMultiplicity(op_ics, xcorners, q0ceed, dm, Qexactloc, Qexact,
165884d34d69SLeila Ghaffari                                restrictq, &ctxSetup, ftime); CHKERRQ(ierr);
1659ccaff030SJeremy L Thompson 
1660ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1661ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1662cfa64770SLeila Ghaffari     CeedVectorDestroy(&q0ceed);
1663ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1664ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1665ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
166684d34d69SLeila Ghaffari     // Clean up vectors
166784d34d69SLeila Ghaffari     ierr = DMRestoreLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
166884d34d69SLeila Ghaffari     ierr = VecDestroy(&Qexact); CHKERRQ(ierr);
1669ccaff030SJeremy L Thompson   }
1670ccaff030SJeremy L Thompson 
1671ccaff030SJeremy L Thompson   // Output Statistics
1672ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts, &steps); CHKERRQ(ierr);
167384d34d69SLeila Ghaffari   if (testChoice == TEST_NONE) {
1674ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1675ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1676ccaff030SJeremy L Thompson                        steps, (double)ftime); CHKERRQ(ierr);
1677ccaff030SJeremy L Thompson   }
167884d34d69SLeila Ghaffari   // Output numerical values from command line
167984d34d69SLeila Ghaffari   ierr = VecViewFromOptions(Q, NULL, "-vec_view"); CHKERRQ(ierr);
168084d34d69SLeila Ghaffari 
168184d34d69SLeila Ghaffari   // compare reference solution values with current run
168284d34d69SLeila Ghaffari   if (testChoice != TEST_NONE) {
168384d34d69SLeila Ghaffari     PetscViewer viewer;
168484d34d69SLeila Ghaffari     // Read reference file
168584d34d69SLeila Ghaffari     Vec Qref;
168684d34d69SLeila Ghaffari     PetscReal error, Qrefnorm;
168784d34d69SLeila Ghaffari     ierr = VecDuplicate(Q, &Qref); CHKERRQ(ierr);
168884d34d69SLeila Ghaffari     ierr = PetscViewerBinaryOpen(comm, test->filepath, FILE_MODE_READ, &viewer);
168984d34d69SLeila Ghaffari     CHKERRQ(ierr);
169084d34d69SLeila Ghaffari     ierr = VecLoad(Qref, viewer); CHKERRQ(ierr);
169184d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
169284d34d69SLeila Ghaffari 
169384d34d69SLeila Ghaffari     // Compute error with respect to reference solution
169484d34d69SLeila Ghaffari     ierr = VecAXPY(Q, -1.0, Qref);  CHKERRQ(ierr);
169584d34d69SLeila Ghaffari     ierr = VecNorm(Qref, NORM_MAX, &Qrefnorm); CHKERRQ(ierr);
169684d34d69SLeila Ghaffari     ierr = VecScale(Q, 1./Qrefnorm); CHKERRQ(ierr);
169784d34d69SLeila Ghaffari     ierr = VecNorm(Q, NORM_MAX, &error); CHKERRQ(ierr);
169884d34d69SLeila Ghaffari     ierr = VecDestroy(&Qref); CHKERRQ(ierr);
169984d34d69SLeila Ghaffari     // Check error
170084d34d69SLeila Ghaffari     if (error > test->testtol) {
170184d34d69SLeila Ghaffari       ierr = PetscPrintf(PETSC_COMM_WORLD,
170284d34d69SLeila Ghaffari                          "Test failed with error norm %g\n",
170384d34d69SLeila Ghaffari                          (double)error); CHKERRQ(ierr);
170484d34d69SLeila Ghaffari     }
170584d34d69SLeila Ghaffari     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
170684d34d69SLeila Ghaffari   }
17079cf88b28Svaleriabarra 
1708ccaff030SJeremy L Thompson   // Clean up libCEED
17098b982baeSLeila Ghaffari   CeedVectorDestroy(&qdata);
1710ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1711ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1712ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1713ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
171484d34d69SLeila Ghaffari   CeedBasisDestroy(&basisq);
171584d34d69SLeila Ghaffari   CeedBasisDestroy(&basisx);
171684d34d69SLeila Ghaffari   CeedBasisDestroy(&basisxc);
171784d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictq);
171884d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictx);
171984d34d69SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdi);
1720ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1721ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1722ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1723ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1724ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1725ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
17266a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
17276a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
17286a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
17296a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
17306a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
17316a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
17326a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
17337ed8b9c7SLeila Ghaffari   CeedQFunctionDestroy(&qf_applySur);
1734ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1735ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1736ccaff030SJeremy L Thompson 
1737ccaff030SJeremy L Thompson   // Clean up PETSc
1738ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1739ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1740ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1741ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1742ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1743ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1744ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1745ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1746ccaff030SJeremy L Thompson   return PetscFinalize();
1747ccaff030SJeremy L Thompson }
1748ccaff030SJeremy L Thompson 
1749