xref: /libCEED/examples/fluids/navierstokes.c (revision 1bb4263918babe26a363fff2527ef14c0f7af2a5)
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 //
31ea6e0f84SLeila Ghaffari //     ./navierstokes -ceed /cpu/self -problem density_current -degree_volume 1
32ea6e0f84SLeila Ghaffari //     ./navierstokes -ceed /gpu/occa -problem advection -degree_volume 1
33ccaff030SJeremy L Thompson //
34ea6e0f84SLeila Ghaffari //TESTARGS -ceed {ceed_resource} -test -degree_volume 1
35ccaff030SJeremy L Thompson 
36ccaff030SJeremy L Thompson /// @file
37ccaff030SJeremy L Thompson /// Navier-Stokes example using PETSc
38ccaff030SJeremy L Thompson 
39ccaff030SJeremy L Thompson const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n";
40ccaff030SJeremy L Thompson 
41ccaff030SJeremy L Thompson #include <petscts.h>
42ccaff030SJeremy L Thompson #include <petscdmplex.h>
43ccaff030SJeremy L Thompson #include <ceed.h>
44ccaff030SJeremy L Thompson #include <stdbool.h>
45ccaff030SJeremy L Thompson #include <petscsys.h>
46ccaff030SJeremy L Thompson #include "common.h"
47b0137797SLeila Ghaffari #include "setup-boundary.h"
48ccaff030SJeremy L Thompson #include "advection.h"
49ccaff030SJeremy L Thompson #include "advection2d.h"
50ccaff030SJeremy L Thompson #include "densitycurrent.h"
51ccaff030SJeremy L Thompson 
52ccaff030SJeremy L Thompson // Problem Options
53ccaff030SJeremy L Thompson typedef enum {
54ccaff030SJeremy L Thompson   NS_DENSITY_CURRENT = 0,
55ccaff030SJeremy L Thompson   NS_ADVECTION = 1,
56ccaff030SJeremy L Thompson   NS_ADVECTION2D = 2,
57ccaff030SJeremy L Thompson } problemType;
58ccaff030SJeremy L Thompson static const char *const problemTypes[] = {
59ccaff030SJeremy L Thompson   "density_current",
60ccaff030SJeremy L Thompson   "advection",
61ccaff030SJeremy L Thompson   "advection2d",
620c6c0b13SLeila Ghaffari   "problemType","NS_",0
63ccaff030SJeremy L Thompson };
64ccaff030SJeremy L Thompson 
65ccaff030SJeremy L Thompson typedef enum {
66ccaff030SJeremy L Thompson   STAB_NONE = 0,
67ccaff030SJeremy L Thompson   STAB_SU = 1,   // Streamline Upwind
68ccaff030SJeremy L Thompson   STAB_SUPG = 2, // Streamline Upwind Petrov-Galerkin
69ccaff030SJeremy L Thompson } StabilizationType;
70ccaff030SJeremy L Thompson static const char *const StabilizationTypes[] = {
710c6c0b13SLeila Ghaffari   "NONE",
72ccaff030SJeremy L Thompson   "SU",
73ccaff030SJeremy L Thompson   "SUPG",
74ccaff030SJeremy L Thompson   "StabilizationType", "STAB_", NULL
75ccaff030SJeremy L Thompson };
76ccaff030SJeremy L Thompson 
77ccaff030SJeremy L Thompson // Problem specific data
78ccaff030SJeremy L Thompson typedef struct {
79*1bb42639SLeila Ghaffari   CeedInt dim, qdatasizeVol;
80ea6e0f84SLeila Ghaffari   CeedQFunctionUser setupVol, setupSur, ics, applyVol_rhs, applySur_rhs,
81ea6e0f84SLeila Ghaffari                     applyVol_ifunction, applySur_ifunction;
82ccaff030SJeremy L Thompson   PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt,
83ccaff030SJeremy L Thompson                        PetscScalar[], void *);
84ea6e0f84SLeila Ghaffari   const char *setupVol_loc, *setupSur_loc, *ics_loc, *applyVol_rhs_loc,
85ea6e0f84SLeila Ghaffari              *applySur_rhs_loc, *applyVol_ifunction_loc, *applySur_ifunction_loc;
86ccaff030SJeremy L Thompson   const bool non_zero_time;
87ccaff030SJeremy L Thompson } problemData;
88ccaff030SJeremy L Thompson 
89ccaff030SJeremy L Thompson problemData problemOptions[] = {
90ccaff030SJeremy L Thompson   [NS_DENSITY_CURRENT] = {
91ccaff030SJeremy L Thompson     .dim                       = 3,
92ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
93b0137797SLeila Ghaffari     .setupVol                  = Setup,
94b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
95356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
96356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
97ccaff030SJeremy L Thompson     .ics                       = ICsDC,
98ccaff030SJeremy L Thompson     .ics_loc                   = ICsDC_loc,
99c96c872fSLeila Ghaffari     .applyVol_rhs              = DC,
100c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = DC_loc,
101ea6e0f84SLeila Ghaffari   //.applySur_rhs              = DC_Sur,
102ea6e0f84SLeila Ghaffari   //.applySur_rhs_loc          = DC_Sur_loc,
103c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_DC,
104c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_DC_loc,
105ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_DC_Sur,
106ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_DC_Sur_loc,
107ccaff030SJeremy L Thompson     .bc                        = Exact_DC,
1080c6c0b13SLeila Ghaffari     .non_zero_time             = false,
109ccaff030SJeremy L Thompson   },
110ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
111ccaff030SJeremy L Thompson     .dim                       = 3,
112ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
113b0137797SLeila Ghaffari     .setupVol                  = Setup,
114b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
115356fbf4bSLeila Ghaffari     .setupSur                  = SetupBoundary,
116356fbf4bSLeila Ghaffari     .setupSur_loc              = SetupBoundary_loc,
117ccaff030SJeremy L Thompson     .ics                       = ICsAdvection,
118ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection_loc,
119c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection,
120c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection_loc,
121ea6e0f84SLeila Ghaffari   //.applySur_rhs              = Advection_Sur,
122ea6e0f84SLeila Ghaffari   //.applySur_rhs_loc          = Advection_Sur_loc,
123c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection,
124c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection_loc,
125ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection_Sur,
126ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection_Sur_loc,
127ccaff030SJeremy L Thompson     .bc                        = Exact_Advection,
1280c6c0b13SLeila Ghaffari     .non_zero_time             = true,
129ccaff030SJeremy L Thompson   },
130ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
131ccaff030SJeremy L Thompson     .dim                       = 2,
132ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 5,
133c96c872fSLeila Ghaffari     .setupVol                  = Setup2d,
134c96c872fSLeila Ghaffari     .setupVol_loc              = Setup2d_loc,
135b0137797SLeila Ghaffari     .setupSur                  = SetupBoundary2d,
136b0137797SLeila Ghaffari     .setupSur_loc              = SetupBoundary2d_loc,
137ccaff030SJeremy L Thompson     .ics                       = ICsAdvection2d,
138ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection2d_loc,
139c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection2d,
140c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection2d_loc,
141b0137797SLeila Ghaffari     .applySur_rhs              = Advection2d_Sur,
142b0137797SLeila Ghaffari     .applySur_rhs_loc          = Advection2d_Sur_loc,
143c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection2d,
144c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection2d_loc,
145ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection2d_Sur,
146ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection2d_Sur_loc,
147ccaff030SJeremy L Thompson     .bc                        = Exact_Advection2d,
1480c6c0b13SLeila Ghaffari     .non_zero_time             = true,
149ccaff030SJeremy L Thompson   },
150ccaff030SJeremy L Thompson };
151ccaff030SJeremy L Thompson 
152ccaff030SJeremy L Thompson // PETSc user data
153ccaff030SJeremy L Thompson typedef struct User_ *User;
154ccaff030SJeremy L Thompson typedef struct Units_ *Units;
155ccaff030SJeremy L Thompson 
156ccaff030SJeremy L Thompson struct User_ {
157ccaff030SJeremy L Thompson   MPI_Comm comm;
158ccaff030SJeremy L Thompson   PetscInt outputfreq;
159ccaff030SJeremy L Thompson   DM dm;
160ccaff030SJeremy L Thompson   DM dmviz;
161ccaff030SJeremy L Thompson   Mat interpviz;
162ccaff030SJeremy L Thompson   Ceed ceed;
163ccaff030SJeremy L Thompson   Units units;
164ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
165ea6e0f84SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs_sur, op_rhs,
166ea6e0f84SLeila Ghaffari                op_ifunction_vol, op_ifunction_sur, op_ifunction;
167ccaff030SJeremy L Thompson   Vec M;
168ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
169ccaff030SJeremy L Thompson   PetscInt contsteps;
170ccaff030SJeremy L Thompson };
171ccaff030SJeremy L Thompson 
172ccaff030SJeremy L Thompson struct Units_ {
173ccaff030SJeremy L Thompson   // fundamental units
174ccaff030SJeremy L Thompson   PetscScalar meter;
175ccaff030SJeremy L Thompson   PetscScalar kilogram;
176ccaff030SJeremy L Thompson   PetscScalar second;
177ccaff030SJeremy L Thompson   PetscScalar Kelvin;
178ccaff030SJeremy L Thompson   // derived units
179ccaff030SJeremy L Thompson   PetscScalar Pascal;
180ccaff030SJeremy L Thompson   PetscScalar JperkgK;
181ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
182ccaff030SJeremy L Thompson   PetscScalar WpermK;
183ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
184ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
185ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
186ccaff030SJeremy L Thompson };
187ccaff030SJeremy L Thompson 
188ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
189ccaff030SJeremy L Thompson struct SimpleBC_ {
190ccaff030SJeremy L Thompson   PetscInt nwall, nslip[3];
1910c6c0b13SLeila Ghaffari   PetscInt walls[10], slips[3][10];
192ccaff030SJeremy L Thompson };
193ccaff030SJeremy L Thompson 
194ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
195ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
196ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
197ccaff030SJeremy L Thompson }
198ccaff030SJeremy L Thompson 
199ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
200ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
2010c6c0b13SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value,
202ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
203ccaff030SJeremy L Thompson 
204ccaff030SJeremy L Thompson   PetscSection   section;
2050c6c0b13SLeila Ghaffari   PetscInt       p, Nelem, Ndof, *erestrict, eoffset, nfields, dim,
2060c6c0b13SLeila Ghaffari                  depth;
2070c6c0b13SLeila Ghaffari   DMLabel depthLabel;
2080c6c0b13SLeila Ghaffari   IS depthIS, iterIS;
2090c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
210ccaff030SJeremy L Thompson   PetscErrorCode ierr;
211ccaff030SJeremy L Thompson   Vec Uloc;
212ccaff030SJeremy L Thompson 
213ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
214ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
215ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm,&section); CHKERRQ(ierr);
216ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
217ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
218ccaff030SJeremy L Thompson   fieldoff[0] = 0;
219ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
220ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
221ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
222ccaff030SJeremy L Thompson   }
223ccaff030SJeremy L Thompson 
2240c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2250c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2260c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2270c6c0b13SLeila Ghaffari   if (domainLabel) {
2280c6c0b13SLeila Ghaffari     IS domainIS;
2290c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2300c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2310c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2320c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2330c6c0b13SLeila Ghaffari   } else {
2340c6c0b13SLeila Ghaffari     iterIS = depthIS;
2350c6c0b13SLeila Ghaffari   }
2360c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2370c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
238ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
2390c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
2400c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
241ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
2420c6c0b13SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
2430c6c0b13SLeila Ghaffari                                    &indices, NULL); CHKERRQ(ierr);
2440c6c0b13SLeila Ghaffari     if (numindices % fieldoff[nfields])
2450c6c0b13SLeila Ghaffari       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
2460c6c0b13SLeila Ghaffari                "Number of closure indices not compatible with Cell %D", c);
247ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
248ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
249ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
250ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
251ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
252ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
253ccaff030SJeremy L Thompson           if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j])
254ccaff030SJeremy L Thompson               != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j)
255ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
256ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
257ccaff030SJeremy L Thompson                      c, i, f, j);
258ccaff030SJeremy L Thompson         }
259ccaff030SJeremy L Thompson       }
260ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
261ccaff030SJeremy L Thompson       PetscInt loc = Involute(indices[i*ncomp[0]]);
2626f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
263ccaff030SJeremy L Thompson     }
2640c6c0b13SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
2650c6c0b13SLeila Ghaffari                                        &indices, NULL); CHKERRQ(ierr);
266ccaff030SJeremy L Thompson   }
2670c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
2680c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
2690c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
270ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
2710c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
2720c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
2730c6c0b13SLeila Ghaffari 
274ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
275ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
276ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
2776f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
2786f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
2796f55dfd5Svaleriabarra                             Erestrict);
280ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
281ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
282ccaff030SJeremy L Thompson }
283ccaff030SJeremy L Thompson 
284c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
285bd910870SLeila Ghaffari static PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt ncompx, CeedInt dim,
286c96c872fSLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value, CeedInt P, CeedInt Q,
287c96c872fSLeila Ghaffari     CeedInt qdatasize, CeedElemRestriction *restrictq,
288c96c872fSLeila Ghaffari     CeedElemRestriction *restrictx, CeedElemRestriction *restrictqdi) {
289c96c872fSLeila Ghaffari 
290c96c872fSLeila Ghaffari   DM dmcoord;
291c96c872fSLeila Ghaffari   CeedInt localNelem;
292c96c872fSLeila Ghaffari   CeedInt Qdim = CeedIntPow(Q, dim);
293c96c872fSLeila Ghaffari   PetscErrorCode ierr;
294c96c872fSLeila Ghaffari 
295c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
296c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
297c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
298c96c872fSLeila Ghaffari   CHKERRQ(ierr);
299c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
300c96c872fSLeila Ghaffari   CHKERRQ(ierr);
301c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
302c96c872fSLeila Ghaffari   CHKERRQ(ierr);
303c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
304c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
305c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
306c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
307c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
308c96c872fSLeila Ghaffari }
309c96c872fSLeila Ghaffari 
310ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
311ccaff030SJeremy L Thompson   PetscErrorCode ierr;
312ccaff030SJeremy L Thompson   PetscInt m;
313ccaff030SJeremy L Thompson 
314ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
315ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
316ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
317ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
318ccaff030SJeremy L Thompson }
319ccaff030SJeremy L Thompson 
320ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
321ccaff030SJeremy L Thompson   PetscErrorCode ierr;
322ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
323ccaff030SJeremy L Thompson   PetscScalar *a;
324ccaff030SJeremy L Thompson 
325ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
326ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
327ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
328ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
3290c6c0b13SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",mpetsc,mceed);
330ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
331ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
332ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
333ccaff030SJeremy L Thompson }
334ccaff030SJeremy L Thompson 
335ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
336ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
337ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
338ccaff030SJeremy L Thompson   PetscErrorCode ierr;
339ccaff030SJeremy L Thompson   Vec Qbc;
340ccaff030SJeremy L Thompson 
341ccaff030SJeremy L Thompson   PetscFunctionBegin;
342ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
343ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
344ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
345ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
346ccaff030SJeremy L Thompson }
347ccaff030SJeremy L Thompson 
348ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
349ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
350ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
351ccaff030SJeremy L Thompson   PetscErrorCode ierr;
352ccaff030SJeremy L Thompson   User user = *(User *)userData;
353ccaff030SJeremy L Thompson   PetscScalar *q, *g;
354ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
355ccaff030SJeremy L Thompson 
356ccaff030SJeremy L Thompson   // Global-to-local
357ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
358ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
359ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
360ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
361ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
362ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
363ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
364ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
365ccaff030SJeremy L Thompson 
366ccaff030SJeremy L Thompson   // Ceed Vectors
367ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
368ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
369ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
370ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
371ccaff030SJeremy L Thompson 
372ccaff030SJeremy L Thompson   // Apply CEED operator
373ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
374ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
375ccaff030SJeremy L Thompson 
376ccaff030SJeremy L Thompson   // Restore vectors
377ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
378ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
379ccaff030SJeremy L Thompson 
380ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
381ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
382ccaff030SJeremy L Thompson 
383ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
384ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
385ccaff030SJeremy L Thompson   CHKERRQ(ierr);
386ccaff030SJeremy L Thompson 
387ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
388ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
389ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
390ccaff030SJeremy L Thompson }
391ccaff030SJeremy L Thompson 
392ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
393ccaff030SJeremy L Thompson                                    void *userData) {
394ccaff030SJeremy L Thompson   PetscErrorCode ierr;
395ccaff030SJeremy L Thompson   User user = *(User *)userData;
396ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
397ccaff030SJeremy L Thompson   PetscScalar *g;
398ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
399ccaff030SJeremy L Thompson 
400ccaff030SJeremy L Thompson   // Global-to-local
401ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
402ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
403ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
404ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
405ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
406ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
407ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
408ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
409ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
410ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
411ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
412ccaff030SJeremy L Thompson 
413ccaff030SJeremy L Thompson   // Ceed Vectors
414ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
415ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
416ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
417ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
418ccaff030SJeremy L Thompson                      (PetscScalar *)q);
419ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
420ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
421ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
422ccaff030SJeremy L Thompson 
423ccaff030SJeremy L Thompson   // Apply CEED operator
424ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
425ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
426ccaff030SJeremy L Thompson 
427ccaff030SJeremy L Thompson   // Restore vectors
428ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
429ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
430ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
431ccaff030SJeremy L Thompson 
432ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
433ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
434ccaff030SJeremy L Thompson 
435ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
436ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
437ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
438ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
439ccaff030SJeremy L Thompson }
440ccaff030SJeremy L Thompson 
441ccaff030SJeremy L Thompson // User provided TS Monitor
442ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
443ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
444ccaff030SJeremy L Thompson   User user = ctx;
445ccaff030SJeremy L Thompson   Vec Qloc;
446ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
447ccaff030SJeremy L Thompson   PetscViewer viewer;
448ccaff030SJeremy L Thompson   PetscErrorCode ierr;
449ccaff030SJeremy L Thompson 
450ccaff030SJeremy L Thompson   // Set up output
451ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
452ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
453ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
454ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
455ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
456ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
457ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
458ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
459ccaff030SJeremy L Thompson 
460ccaff030SJeremy L Thompson   // Output
461ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
462ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
463ccaff030SJeremy L Thompson   CHKERRQ(ierr);
464ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
465ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
466ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
467ccaff030SJeremy L Thompson   if (user->dmviz) {
468ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
469ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
470ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
471ccaff030SJeremy L Thompson 
472ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
473ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
474ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
475ccaff030SJeremy L Thompson     CHKERRQ(ierr);
476ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
477ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
478ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
479ccaff030SJeremy L Thompson     CHKERRQ(ierr);
480ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
481ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
482ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
483ccaff030SJeremy L Thompson     CHKERRQ(ierr);
484ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
485ccaff030SJeremy L Thompson                               filepath_refined,
486ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
487ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
488ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
489ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
490ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
491ccaff030SJeremy L Thompson   }
4920c6c0b13SLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
493ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
494ccaff030SJeremy L Thompson 
495ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
496ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
497ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
499ccaff030SJeremy L Thompson   CHKERRQ(ierr);
500ccaff030SJeremy L Thompson   ierr = VecView(Q, viewer); CHKERRQ(ierr);
501ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
502ccaff030SJeremy L Thompson 
503ccaff030SJeremy L Thompson   // Save time stamp
504ccaff030SJeremy L Thompson   // Dimensionalize time back
505ccaff030SJeremy L Thompson   time /= user->units->second;
506ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
507ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
508ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
509ccaff030SJeremy L Thompson   CHKERRQ(ierr);
510ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
511ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
512ccaff030SJeremy L Thompson   #else
513ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
514ccaff030SJeremy L Thompson   #endif
515ccaff030SJeremy L Thompson   CHKERRQ(ierr);
516ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
517ccaff030SJeremy L Thompson 
518ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
519ccaff030SJeremy L Thompson }
520ccaff030SJeremy L Thompson 
5210c6c0b13SLeila Ghaffari static PetscErrorCode ICs_PetscMultiplicity(CeedOperator op_ics,
522ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
523ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
524ccaff030SJeremy L Thompson   PetscErrorCode ierr;
525ccaff030SJeremy L Thompson   CeedVector multlvec;
526ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
527ccaff030SJeremy L Thompson 
528ccaff030SJeremy L Thompson   ctxSetup->time = time;
529ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
530ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
531ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
532ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
533ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
534ccaff030SJeremy L Thompson 
535ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
536ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
537ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
538ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
539ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
540ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
541ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
542ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
543ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
544ccaff030SJeremy L Thompson   CHKERRQ(ierr);
545ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
546ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
547ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
548ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
549ccaff030SJeremy L Thompson 
550ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
551ccaff030SJeremy L Thompson }
552ccaff030SJeremy L Thompson 
553ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
554ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
555ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
556ccaff030SJeremy L Thompson   PetscErrorCode ierr;
557ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
558ccaff030SJeremy L Thompson   CeedOperator op_mass;
559ccaff030SJeremy L Thompson   CeedVector mceed;
560ccaff030SJeremy L Thompson   Vec Mloc;
561ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
562ccaff030SJeremy L Thompson 
563ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
564ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
565ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
566ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
567ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
568ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
569ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
570ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
571ccaff030SJeremy L Thompson 
572ccaff030SJeremy L Thompson   // Create the mass operator
573ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
574ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
575ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
576ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
577ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
578ccaff030SJeremy L Thompson 
579ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
580ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
581ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
582ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
583ccaff030SJeremy L Thompson 
584ccaff030SJeremy L Thompson   {
585ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
586ccaff030SJeremy L Thompson     CeedVector onesvec;
587ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
588ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
589ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
590ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
591ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
592ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
593ccaff030SJeremy L Thompson   }
594ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
595ccaff030SJeremy L Thompson 
596ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
597ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
598ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
599ccaff030SJeremy L Thompson 
600ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
601ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
602ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
603ccaff030SJeremy L Thompson }
604ccaff030SJeremy L Thompson 
6050c6c0b13SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
606ff6701fcSJed Brown                        SimpleBC bc, void *ctxSetup) {
607ccaff030SJeremy L Thompson   PetscErrorCode ierr;
608ccaff030SJeremy L Thompson 
609ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
610ccaff030SJeremy L Thompson   {
611ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
612ccaff030SJeremy L Thompson     PetscFE fe;
6130c6c0b13SLeila Ghaffari     PetscSpace fespace;
614ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
615ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
616ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
61732ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
618ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
619ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
620ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
6210c6c0b13SLeila Ghaffari     /* Wall boundary conditions are zero velocity and zero flux for density and energy */
6220c6c0b13SLeila Ghaffari     {
6230c6c0b13SLeila Ghaffari       PetscInt comps[3] = {1, 2, 3};
6240c6c0b13SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
6250c6c0b13SLeila Ghaffari                            3, comps, (void(*)(void))problem->bc,
6260c6c0b13SLeila Ghaffari                            bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
6270c6c0b13SLeila Ghaffari     }
62807af6069Svaleriabarra     {
62907af6069Svaleriabarra       PetscInt comps[1] = {1};
63007af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
63107af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
63207af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
63307af6069Svaleriabarra       comps[0] = 2;
63407af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
63507af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
63607af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
63707af6069Svaleriabarra       comps[0] = 3;
63807af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
63907af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
64007af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
64107af6069Svaleriabarra     }
642ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE,NULL);
643ccaff030SJeremy L Thompson     CHKERRQ(ierr);
6440c6c0b13SLeila Ghaffari     ierr = PetscFEGetBasisSpace(fe, &fespace); CHKERRQ(ierr);
645ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
646ccaff030SJeremy L Thompson   }
647ccaff030SJeremy L Thompson   {
648ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
649ccaff030SJeremy L Thompson     PetscSection section;
650ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
651ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
652ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
653ccaff030SJeremy L Thompson     CHKERRQ(ierr);
654ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
655ccaff030SJeremy L Thompson     CHKERRQ(ierr);
656ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
657ccaff030SJeremy L Thompson     CHKERRQ(ierr);
658ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
659ccaff030SJeremy L Thompson     CHKERRQ(ierr);
660ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
661ccaff030SJeremy L Thompson     CHKERRQ(ierr);
662ccaff030SJeremy L Thompson   }
663ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
664ccaff030SJeremy L Thompson }
665ccaff030SJeremy L Thompson 
666ccaff030SJeremy L Thompson int main(int argc, char **argv) {
667ccaff030SJeremy L Thompson   PetscInt ierr;
668ccaff030SJeremy L Thompson   MPI_Comm comm;
6690c6c0b13SLeila Ghaffari   DM dm, dmcoord, dmviz, dmvizfine;
670ccaff030SJeremy L Thompson   Mat interpviz;
671ccaff030SJeremy L Thompson   TS ts;
672ccaff030SJeremy L Thompson   TSAdapt adapt;
673ccaff030SJeremy L Thompson   User user;
674ccaff030SJeremy L Thompson   Units units;
675ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
676c96c872fSLeila Ghaffari   PetscInt localNelemVol, localNelemSur, lnodes, steps;
677ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
678ccaff030SJeremy L Thompson   PetscMPIInt rank;
679ccaff030SJeremy L Thompson   PetscScalar ftime;
680ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
681ccaff030SJeremy L Thompson   Ceed ceed;
682ea6e0f84SLeila Ghaffari   CeedInt numP_Vol, numP_Sur, numQ_Vol, numQ_Sur;
683c96c872fSLeila Ghaffari   CeedVector xcorners, qdataVol, qdataSur, q0ceedVol, q0ceedSur;
684ea6e0f84SLeila Ghaffari   CeedBasis basisxVol, basisxcVol, basisqVol, basisxSur, basisxcSur, basisqSur;
685bd910870SLeila Ghaffari   CeedElemRestriction restrictxVol, restrictqVol, restrictqdiVol,
686bd910870SLeila Ghaffari                       restrictxSur, restrictqSur, restrictqdiSur;
687ea6e0f84SLeila Ghaffari   CeedQFunction qf_setupVol, qf_setupSur, qf_ics, qf_rhsVol, qf_rhsSur,
688ea6e0f84SLeila Ghaffari                 qf_ifunctionVol, qf_ifunctionSur;
689ea6e0f84SLeila Ghaffari   CeedOperator op_setupVol, op_setupSur, op_ics;
690ccaff030SJeremy L Thompson   CeedScalar Rd;
691ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
692ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
693ccaff030SJeremy L Thompson   problemType problemChoice;
694ccaff030SJeremy L Thompson   problemData *problem = NULL;
695ccaff030SJeremy L Thompson   StabilizationType stab;
6960c6c0b13SLeila Ghaffari   PetscBool   test, implicit;
697cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
698ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
6990c6c0b13SLeila Ghaffari     .nwall = 6,
7000c6c0b13SLeila Ghaffari     .walls = {1,2,3,4,5,6},
701ccaff030SJeremy L Thompson   };
702ccaff030SJeremy L Thompson   double start, cpu_time_used;
703ccaff030SJeremy L Thompson 
704ccaff030SJeremy L Thompson   // Create the libCEED contexts
705ccaff030SJeremy L Thompson   PetscScalar meter     = 1e-2;     // 1 meter in scaled length units
706ccaff030SJeremy L Thompson   PetscScalar second    = 1e-2;     // 1 second in scaled time units
707ccaff030SJeremy L Thompson   PetscScalar kilogram  = 1e-6;     // 1 kilogram in scaled mass units
708ccaff030SJeremy L Thompson   PetscScalar Kelvin    = 1;        // 1 Kelvin in scaled temperature units
709ccaff030SJeremy L Thompson   CeedScalar theta0     = 300.;     // K
710ccaff030SJeremy L Thompson   CeedScalar thetaC     = -15.;     // K
711ccaff030SJeremy L Thompson   CeedScalar P0         = 1.e5;     // Pa
712ccaff030SJeremy L Thompson   CeedScalar N          = 0.01;     // 1/s
713ccaff030SJeremy L Thompson   CeedScalar cv         = 717.;     // J/(kg K)
714ccaff030SJeremy L Thompson   CeedScalar cp         = 1004.;    // J/(kg K)
715ccaff030SJeremy L Thompson   CeedScalar g          = 9.81;     // m/s^2
716ccaff030SJeremy L Thompson   CeedScalar lambda     = -2./3.;   // -
717ccaff030SJeremy L Thompson   CeedScalar mu         = 75.;      // Pa s, dynamic viscosity
718ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
719ccaff030SJeremy L Thompson   CeedScalar k          = 0.02638;  // W/(m K)
720ccaff030SJeremy L Thompson   CeedScalar CtauS      = 0.;       // dimensionless
721ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
722ccaff030SJeremy L Thompson   PetscScalar lx        = 8000.;    // m
723ccaff030SJeremy L Thompson   PetscScalar ly        = 8000.;    // m
724ccaff030SJeremy L Thompson   PetscScalar lz        = 4000.;    // m
725ccaff030SJeremy L Thompson   CeedScalar rc         = 1000.;    // m (Radius of bubble)
726ccaff030SJeremy L Thompson   PetscScalar resx      = 1000.;    // m (resolution in x)
727ccaff030SJeremy L Thompson   PetscScalar resy      = 1000.;    // m (resolution in y)
728ccaff030SJeremy L Thompson   PetscScalar resz      = 1000.;    // m (resolution in z)
729ccaff030SJeremy L Thompson   PetscInt outputfreq   = 10;       // -
730ccaff030SJeremy L Thompson   PetscInt contsteps    = 0;        // -
731ea6e0f84SLeila Ghaffari   PetscInt degreeVol    = 1;        // -
732ea6e0f84SLeila Ghaffari   PetscInt degreeSur    = 1;        // -
733ea6e0f84SLeila Ghaffari   PetscInt qextraVol    = 2;        // -
734ea6e0f84SLeila Ghaffari   PetscInt qextraSur    = 2;        // -
7350c6c0b13SLeila Ghaffari   DMBoundaryType periodicity[] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
7360c6c0b13SLeila Ghaffari                                   DM_BOUNDARY_NONE
7370c6c0b13SLeila Ghaffari                                  };
738ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
739ccaff030SJeremy L Thompson 
740ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
741ccaff030SJeremy L Thompson   if (ierr) return ierr;
742ccaff030SJeremy L Thompson 
743ccaff030SJeremy L Thompson   // Allocate PETSc context
744ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
745ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
746ccaff030SJeremy L Thompson 
747ccaff030SJeremy L Thompson   // Parse command line options
748ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
749ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
750ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
751ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
752ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
753ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
7540c6c0b13SLeila Ghaffari   ierr = PetscOptionsBool("-test", "Run in test mode",
7550c6c0b13SLeila Ghaffari                           NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr);
756ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
757ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
758ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
759ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
760ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
761ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
762ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
763ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
764ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
765ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
766ccaff030SJeremy L Thompson   CHKERRQ(ierr);
767ccaff030SJeremy L Thompson   {
768ccaff030SJeremy L Thompson     PetscInt len;
769ccaff030SJeremy L Thompson     PetscBool flg;
770ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
771ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
772ccaff030SJeremy L Thompson                                 NULL, bc.walls,
773ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
774ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
775ccaff030SJeremy L Thompson     if (flg) bc.nwall = len;
776ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
777ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
778ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
779ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
780ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
781ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
782ccaff030SJeremy L Thompson                                    &len), &flg);
783ccaff030SJeremy L Thompson       CHKERRQ(ierr);
7840c6c0b13SLeila Ghaffari       if (flg) bc.nslip[j] = len;
785ccaff030SJeremy L Thompson     }
786ccaff030SJeremy L Thompson   }
787cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
788cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
789cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
790ccaff030SJeremy L Thompson   CHKERRQ(ierr);
791ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
792ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
793ccaff030SJeremy L Thompson   meter = fabs(meter);
794ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
795ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
796ccaff030SJeremy L Thompson   second = fabs(second);
797ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
798ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
799ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
800ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
801ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
802ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
803ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
804ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
805ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
806ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
807ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
808ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
809ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
810ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
811ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
812ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
813ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
814ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
815ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
816ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
817ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
818ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
819ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
820ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
821ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
822ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
823ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
824ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
825ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
826ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
827ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
828ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
829ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
830ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
831ccaff030SJeremy L Thompson   CHKERRQ(ierr);
832ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
833ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
834ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
835ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
836ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
837ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
838ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
839ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
840ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
841ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
842ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
843ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
844ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
845ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
846ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
8470c6c0b13SLeila Ghaffari   ierr = PetscOptionsEnumArray("-periodicity", "Periodicity per direction",
8480c6c0b13SLeila Ghaffari                                NULL, DMBoundaryTypes, (PetscEnum *)periodicity,
8490c6c0b13SLeila Ghaffari                                &n, NULL); CHKERRQ(ierr);
8500c6c0b13SLeila Ghaffari   n = problem->dim;
851ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
852ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
853ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
854ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
855ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
856ccaff030SJeremy L Thompson   n = problem->dim;
857ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
858ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
859ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
860ccaff030SJeremy L Thompson   {
861ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
862ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
863ccaff030SJeremy L Thompson     if (norm > 0) {
864ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
865ccaff030SJeremy L Thompson     }
866ccaff030SJeremy L Thompson   }
867ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
868ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
869ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
870ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
871ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
872ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-degree_volume", "Polynomial degree of finite elements for the Volume",
873ea6e0f84SLeila Ghaffari                          NULL, degreeVol, &degreeVol, NULL); CHKERRQ(ierr);
874ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-qextra_volume", "Number of extra quadrature points for the Volume",
875ea6e0f84SLeila Ghaffari                          NULL, qextraVol, &qextraVol, NULL); CHKERRQ(ierr);
8760c6c0b13SLeila Ghaffari   PetscStrncpy(user->outputfolder, ".", 2);
877ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
878ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
879ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
880ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
881ccaff030SJeremy L Thompson 
882ccaff030SJeremy L Thompson   // Define derived units
883ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
884ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
885ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
886ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
887ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
888ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
889ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
890ccaff030SJeremy L Thompson 
891ccaff030SJeremy L Thompson   // Scale variables to desired units
892ccaff030SJeremy L Thompson   theta0 *= Kelvin;
893ccaff030SJeremy L Thompson   thetaC *= Kelvin;
894ccaff030SJeremy L Thompson   P0 *= Pascal;
895ccaff030SJeremy L Thompson   N *= (1./second);
896ccaff030SJeremy L Thompson   cv *= JperkgK;
897ccaff030SJeremy L Thompson   cp *= JperkgK;
898ccaff030SJeremy L Thompson   Rd = cp - cv;
899ccaff030SJeremy L Thompson   g *= mpersquareds;
900ccaff030SJeremy L Thompson   mu *= Pascal * second;
901ccaff030SJeremy L Thompson   k *= WpermK;
902ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
903ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
904ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
905ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
906ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
907ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
908ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
909ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
910ccaff030SJeremy L Thompson 
911ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
912*1bb42639SLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol;
913ccaff030SJeremy L Thompson   // Set up the libCEED context
914ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup =  {
915ccaff030SJeremy L Thompson     .theta0 = theta0,
916ccaff030SJeremy L Thompson     .thetaC = thetaC,
917ccaff030SJeremy L Thompson     .P0 = P0,
918ccaff030SJeremy L Thompson     .N = N,
919ccaff030SJeremy L Thompson     .cv = cv,
920ccaff030SJeremy L Thompson     .cp = cp,
921ccaff030SJeremy L Thompson     .Rd = Rd,
922ccaff030SJeremy L Thompson     .g = g,
923ccaff030SJeremy L Thompson     .rc = rc,
924ccaff030SJeremy L Thompson     .lx = lx,
925ccaff030SJeremy L Thompson     .ly = ly,
926ccaff030SJeremy L Thompson     .lz = lz,
9270c6c0b13SLeila Ghaffari     .periodicity0 = periodicity[0],
9280c6c0b13SLeila Ghaffari     .periodicity1 = periodicity[1],
9290c6c0b13SLeila Ghaffari     .periodicity2 = periodicity[2],
930ccaff030SJeremy L Thompson     .center[0] = center[0],
931ccaff030SJeremy L Thompson     .center[1] = center[1],
932ccaff030SJeremy L Thompson     .center[2] = center[2],
933ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
934ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
935ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
936ccaff030SJeremy L Thompson     .time = 0,
937ccaff030SJeremy L Thompson   };
938ccaff030SJeremy L Thompson 
939ccaff030SJeremy L Thompson   {
940ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
941ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
9420c6c0b13SLeila Ghaffari                                periodicity, PETSC_TRUE, &dm);
943ccaff030SJeremy L Thompson     CHKERRQ(ierr);
944ccaff030SJeremy L Thompson   }
9450c6c0b13SLeila Ghaffari   if (1) {
946ccaff030SJeremy L Thompson     DM               dmDist = NULL;
947ccaff030SJeremy L Thompson     PetscPartitioner part;
948ccaff030SJeremy L Thompson 
949ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
950ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
951ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
952ccaff030SJeremy L Thompson     if (dmDist) {
953ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
954ccaff030SJeremy L Thompson       dm  = dmDist;
955ccaff030SJeremy L Thompson     }
956ccaff030SJeremy L Thompson   }
957ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
958ccaff030SJeremy L Thompson 
959ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
960ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
961ea6e0f84SLeila Ghaffari   ierr = SetUpDM(dm, problem, degreeVol, &bc, &ctxSetup); CHKERRQ(ierr);
9620c6c0b13SLeila Ghaffari   if (!test) {
9630c6c0b13SLeila Ghaffari     ierr = PetscPrintf(PETSC_COMM_WORLD,
964ea6e0f84SLeila Ghaffari                        "Degree of Volumetric FEM Space: %D\n",
965ea6e0f84SLeila Ghaffari                        degreeVol); CHKERRQ(ierr);
9660c6c0b13SLeila Ghaffari   }
967ccaff030SJeremy L Thompson   dmviz = NULL;
968ccaff030SJeremy L Thompson   interpviz = NULL;
969ccaff030SJeremy L Thompson   if (viz_refine) {
970ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
971ff6701fcSJed Brown 
972ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
973ff6701fcSJed Brown     dmhierarchy[0] = dm;
974ea6e0f84SLeila Ghaffari     for (PetscInt i = 0, d = degreeVol; i < viz_refine; i++) {
975ff6701fcSJed Brown       Mat interp_next;
976ff6701fcSJed Brown 
977ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
978ccaff030SJeremy L Thompson       CHKERRQ(ierr);
979ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
980ff6701fcSJed Brown       d = (d + 1) / 2;
981ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
982ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
983ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
984ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
985ff6701fcSJed Brown       if (!i) interpviz = interp_next;
986ff6701fcSJed Brown       else {
987ff6701fcSJed Brown         Mat C;
988ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
989ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
990ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
991ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
992ff6701fcSJed Brown         interpviz = C;
993ff6701fcSJed Brown       }
994ff6701fcSJed Brown     }
995cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
996ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
997cb3e2689Svaleriabarra     }
998ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
999ccaff030SJeremy L Thompson   }
1000ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1001ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1002ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1003ccaff030SJeremy L Thompson   lnodes /= ncompq;
1004ccaff030SJeremy L Thompson 
10050c6c0b13SLeila Ghaffari   {
10060c6c0b13SLeila Ghaffari     // Print grid information
1007ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1008ccaff030SJeremy L Thompson     int comm_size;
1009ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1010ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1011ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
1012ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1013ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1014ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
10150c6c0b13SLeila Ghaffari     if (!test) {
10160c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Global FEM dofs: %D (%D owned) on %d ranks\n",
10170c6c0b13SLeila Ghaffari                          gdofs, odofs, comm_size); CHKERRQ(ierr);
10180c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Local FEM nodes: %D\n", lnodes); CHKERRQ(ierr);
10190c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "dm_plex_box_faces: %s\n", box_faces_str);
1020ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1021ccaff030SJeremy L Thompson     }
1022ccaff030SJeremy L Thompson 
10230c6c0b13SLeila Ghaffari   }
10240c6c0b13SLeila Ghaffari 
1025ccaff030SJeremy L Thompson   // Set up global mass vector
1026ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q,&user->M); CHKERRQ(ierr);
1027ccaff030SJeremy L Thompson 
10280c6c0b13SLeila Ghaffari   // Set up CEED
1029ccaff030SJeremy L Thompson   // CEED Bases
1030ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
1031ea6e0f84SLeila Ghaffari   numP_Vol = degreeVol + 1;
1032ea6e0f84SLeila Ghaffari   numQ_Vol = numP_Vol + qextraVol;
1033ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP_Vol, numQ_Vol, CEED_GAUSS,
1034ea6e0f84SLeila Ghaffari                                   &basisqVol);
1035ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ_Vol, CEED_GAUSS,
1036ea6e0f84SLeila Ghaffari                                   &basisxVol);
1037ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP_Vol,
1038ea6e0f84SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcVol);
1039ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1040ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1041ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1042ccaff030SJeremy L Thompson 
1043ccaff030SJeremy L Thompson   // CEED Restrictions
1044c96c872fSLeila Ghaffari   // Restrictions on the Volume
10456a0edaf9SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, ncompx, dim, 0, 0, 0, numP_Vol, numQ_Vol,
10466a0edaf9SLeila Ghaffari                                  qdatasizeVol, &restrictqVol, &restrictxVol,
10476a0edaf9SLeila Ghaffari                                  &restrictqdiVol); CHKERRQ(ierr);
1048ccaff030SJeremy L Thompson 
1049ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1050ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1051ccaff030SJeremy L Thompson 
1052ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1053bd910870SLeila Ghaffari   CeedInt NqptsVol;
1054c96c872fSLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqVol, &NqptsVol);
1055c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictqVol, &localNelemVol);
1056c96c872fSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdataVol);
1057c96c872fSLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &q0ceedVol, NULL);
1058ccaff030SJeremy L Thompson 
1059ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1060ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1061ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1062ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1063ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
1064ea6e0f84SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1065ccaff030SJeremy L Thompson 
1066ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1067ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1068ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1069ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1070ccaff030SJeremy L Thompson 
1071ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1072ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1073ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1074ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1075ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1076ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1077ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1078ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1079ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1080ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1081ccaff030SJeremy L Thompson   }
1082ccaff030SJeremy L Thompson 
1083ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1084ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1085ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1086ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1087ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1088ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1089ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
1090ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1091ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1092ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1093ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1094ccaff030SJeremy L Thompson   }
1095ccaff030SJeremy L Thompson 
1096ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1097ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
1098ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictxVol, basisxVol, CEED_VECTOR_ACTIVE);
1099ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
1100ea6e0f84SLeila Ghaffari                        basisxVol, CEED_VECTOR_NONE);
1101ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdataVol", restrictqdiVol,
1102ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1103ccaff030SJeremy L Thompson 
1104ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1105ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
1106ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictxVol, basisxcVol, CEED_VECTOR_ACTIVE);
1107ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictqVol,
1108ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1109ccaff030SJeremy L Thompson 
1110ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qceed, NULL);
1111ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qdotceed, NULL);
1112ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->gceed, NULL);
1113ccaff030SJeremy L Thompson 
1114ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1115ccaff030SJeremy L Thompson     CeedOperator op;
1116ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
1117ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1118ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1119ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdataVol", restrictqdiVol,
1120ea6e0f84SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataVol);
1121ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1122ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1123ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1124ccaff030SJeremy L Thompson     user->op_rhs = op;
1125ccaff030SJeremy L Thompson   }
1126ccaff030SJeremy L Thompson 
1127ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1128ccaff030SJeremy L Thompson     CeedOperator op;
1129ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
1130ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1131ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1132ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictqVol, basisqVol, user->qdotceed);
1133ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdataVol", restrictqdiVol,
1134ea6e0f84SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataVol);
1135ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1136ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1137ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1138ccaff030SJeremy L Thompson     user->op_ifunction = op;
1139ccaff030SJeremy L Thompson   }
1140ccaff030SJeremy L Thompson 
11416a0edaf9SLeila Ghaffari   //**************************************************************************************//
11426a0edaf9SLeila Ghaffari   // Add boundary Integral (TODO: create a function for all faces)
11436a0edaf9SLeila Ghaffari   //--------------------------------------------------------------------------------------//
11446a0edaf9SLeila Ghaffari   // Set up CEED for the boundaries
11456a0edaf9SLeila Ghaffari   // CEED bases
11466a0edaf9SLeila Ghaffari   CeedInt height = 1;
11476a0edaf9SLeila Ghaffari   CeedInt dimSur = dim - height;
11486a0edaf9SLeila Ghaffari   numP_Sur = degreeSur + 1;
11496a0edaf9SLeila Ghaffari   numQ_Sur = numP_Sur + qextraSur;
11506a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
11516a0edaf9SLeila Ghaffari                                   &basisqSur);
11526a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numQ_Sur, CEED_GAUSS,
11536a0edaf9SLeila Ghaffari                                   &basisxSur);
11546a0edaf9SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dimSur, ncompx, 2, numP_Sur,
11556a0edaf9SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
11566a0edaf9SLeila Ghaffari   // CEED Restrictions
11576a0edaf9SLeila Ghaffari   // Restriction on one face
11586a0edaf9SLeila Ghaffari   DMLabel domainLabel;
1159*1bb42639SLeila Ghaffari   CeedInt qdatasizeSur = 1;
11606a0edaf9SLeila Ghaffari   ierr = DMGetLabel(dm, "Face Sets", &domainLabel); CHKERRQ(ierr);
11616a0edaf9SLeila Ghaffari   ierr = GetRestrictionForDomain(ceed, dm, ncompx, dimSur, height, domainLabel, 2, numP_Sur,
11626a0edaf9SLeila Ghaffari                                  numQ_Sur, qdatasizeSur, &restrictqSur, &restrictxSur,
11636a0edaf9SLeila Ghaffari                                  &restrictqdiSur); CHKERRQ(ierr);
11646a0edaf9SLeila Ghaffari 
11656a0edaf9SLeila Ghaffari   // Create the CEED vectors that will be needed in setup
11666a0edaf9SLeila Ghaffari   CeedInt NqptsSur;
11676a0edaf9SLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
11686a0edaf9SLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictqSur, &localNelemSur);
11696a0edaf9SLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeSur*localNelemSur*NqptsSur, &qdataSur);
11706a0edaf9SLeila Ghaffari 
11716a0edaf9SLeila Ghaffari   // Create the Q-function that builds the quadrature data for the NS operator
11726a0edaf9SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupSur, problem->setupSur_loc,
11736a0edaf9SLeila Ghaffari                               &qf_setupSur);
11746a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "dx", ncompx*dimSur, CEED_EVAL_GRAD);
11756a0edaf9SLeila Ghaffari   CeedQFunctionAddInput(qf_setupSur, "weight", 1, CEED_EVAL_WEIGHT);
11766a0edaf9SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
11776a0edaf9SLeila Ghaffari 
11786a0edaf9SLeila Ghaffari 
11796a0edaf9SLeila Ghaffari   qf_rhsSur = NULL;
11806a0edaf9SLeila Ghaffari   if (problem->applySur_rhs) { // Create the Q-function that defines the action of the RHS operator
11816a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_rhs,
11826a0edaf9SLeila Ghaffari                                 problem->applySur_rhs_loc, &qf_rhsSur);
11836a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "q", ncompq, CEED_EVAL_INTERP);
11846a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "dq", ncompq*dimSur, CEED_EVAL_GRAD);
11856a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
11866a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsSur, "x", ncompx, CEED_EVAL_INTERP);
11876a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsSur, "v", ncompq, CEED_EVAL_INTERP);
11886a0edaf9SLeila Ghaffari   }
11896a0edaf9SLeila Ghaffari 
11906a0edaf9SLeila Ghaffari   qf_ifunctionSur = NULL;
11916a0edaf9SLeila Ghaffari   if (problem->applySur_ifunction) { // Create the Q-function that defines the action of the IFunction
11926a0edaf9SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applySur_ifunction,
11936a0edaf9SLeila Ghaffari                                 problem->applySur_ifunction_loc, &qf_ifunctionSur);
11946a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "q", ncompq, CEED_EVAL_INTERP);
11956a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "dq", ncompq*dimSur, CEED_EVAL_GRAD);
11966a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "qdataSur", qdatasizeSur, CEED_EVAL_NONE);
11976a0edaf9SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionSur, "x", ncompx, CEED_EVAL_INTERP);
11986a0edaf9SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionSur, "v", ncompq, CEED_EVAL_INTERP);
11996a0edaf9SLeila Ghaffari   }
12006a0edaf9SLeila Ghaffari 
12016a0edaf9SLeila Ghaffari   // Create the operator that builds the quadrature data for the NS operator
12026a0edaf9SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupSur, NULL, NULL, &op_setupSur);
12036a0edaf9SLeila Ghaffari   CeedOperatorSetField(op_setupSur, "dx", restrictxSur, basisxSur, CEED_VECTOR_ACTIVE);
12046a0edaf9SLeila Ghaffari   CeedOperatorSetField(op_setupSur, "weight", CEED_ELEMRESTRICTION_NONE,
12056a0edaf9SLeila Ghaffari                        basisxSur, CEED_VECTOR_NONE);
12066a0edaf9SLeila Ghaffari   CeedOperatorSetField(op_setupSur, "qdataSur", restrictqdiSur,
12076a0edaf9SLeila Ghaffari                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
12086a0edaf9SLeila Ghaffari 
12096a0edaf9SLeila Ghaffari 
12106a0edaf9SLeila Ghaffari   if (qf_rhsSur) { // Create the RHS physics operator
12116a0edaf9SLeila Ghaffari     CeedOperator op;
12126a0edaf9SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsSur, NULL, NULL, &op);
12136a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
12146a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
12156a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "qdataSur", restrictqdiSur,
12166a0edaf9SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataSur);
12176a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxSur, basisxSur, xcorners);
12186a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
12196a0edaf9SLeila Ghaffari     user->op_rhs_sur = op;
12206a0edaf9SLeila Ghaffari   }
12216a0edaf9SLeila Ghaffari 
12226a0edaf9SLeila Ghaffari   if (qf_ifunctionSur) { // Create the IFunction operator
12236a0edaf9SLeila Ghaffari     CeedOperator op;
12246a0edaf9SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionSur, NULL, NULL, &op);
12256a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
12266a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
12276a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "qdataSur", restrictqdiSur,
12286a0edaf9SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataSur);
12296a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxSur, basisxSur, xcorners);
12306a0edaf9SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqSur, basisqSur, CEED_VECTOR_ACTIVE);
12316a0edaf9SLeila Ghaffari     user->op_ifunction_sur = op;
12326a0edaf9SLeila Ghaffari   }
12336a0edaf9SLeila Ghaffari   // Composite Operaters
12346a0edaf9SLeila Ghaffari   if (user->op_ifunction_vol) {
12356a0edaf9SLeila Ghaffari     if (user->op_ifunction_sur) {
12366a0edaf9SLeila Ghaffari       // Composite Operators for the IFunction
12376a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_ifunction);
12386a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_vol);
12396a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_sur);
12406a0edaf9SLeila Ghaffari   } else {
12416a0edaf9SLeila Ghaffari     user->op_ifunction = user->op_ifunction_vol;
12426a0edaf9SLeila Ghaffari     }
12436a0edaf9SLeila Ghaffari   }
12446a0edaf9SLeila Ghaffari   if (user->op_rhs_vol) {
12456a0edaf9SLeila Ghaffari     if (user->op_rhs_sur) {
12466a0edaf9SLeila Ghaffari       // Composite Operators for the RHS
12476a0edaf9SLeila Ghaffari       CeedCompositeOperatorCreate(ceed, &user->op_rhs);
12486a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_vol);
12496a0edaf9SLeila Ghaffari       CeedCompositeOperatorAddSub(user->op_rhs, user->op_rhs_sur);
12506a0edaf9SLeila Ghaffari   } else {
12516a0edaf9SLeila Ghaffari     user->op_rhs = user->op_rhs_vol;
12526a0edaf9SLeila Ghaffari     }
12536a0edaf9SLeila Ghaffari   }
12546a0edaf9SLeila Ghaffari   //**************************************************************************************//
1255ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1256ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1257ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1258ccaff030SJeremy L Thompson     .CtauS = CtauS,
1259ccaff030SJeremy L Thompson     .strong_form = strong_form,
1260ccaff030SJeremy L Thompson     .stabilization = stab,
1261ccaff030SJeremy L Thompson   };
1262ccaff030SJeremy L Thompson   switch (problemChoice) {
1263ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1264ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1265ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1266ccaff030SJeremy L Thompson           sizeof ctxNS);
12676a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxNS, sizeof ctxNS);
12686a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxNS,
12696a0edaf9SLeila Ghaffari           sizeof ctxNS);
1270ccaff030SJeremy L Thompson     break;
1271ccaff030SJeremy L Thompson   case NS_ADVECTION:
1272ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1273ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1274ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1275ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1276ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
12776a0edaf9SLeila Ghaffari     if (qf_rhsSur) CeedQFunctionSetContext(qf_rhsSur, &ctxAdvection2d,
12786a0edaf9SLeila Ghaffari                                           sizeof ctxAdvection2d);
12796a0edaf9SLeila Ghaffari     if (qf_ifunctionSur) CeedQFunctionSetContext(qf_ifunctionSur, &ctxAdvection2d,
12806a0edaf9SLeila Ghaffari           sizeof ctxAdvection2d);
1281ccaff030SJeremy L Thompson   }
1282ccaff030SJeremy L Thompson 
1283ccaff030SJeremy L Thompson   // Set up PETSc context
1284ccaff030SJeremy L Thompson   // Set up units structure
1285ccaff030SJeremy L Thompson   units->meter = meter;
1286ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1287ccaff030SJeremy L Thompson   units->second = second;
1288ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1289ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1290ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1291ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1292ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1293ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1294ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1295ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1296ccaff030SJeremy L Thompson 
1297ccaff030SJeremy L Thompson   // Set up user structure
1298ccaff030SJeremy L Thompson   user->comm = comm;
1299ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1300ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1301ccaff030SJeremy L Thompson   user->units = units;
1302ccaff030SJeremy L Thompson   user->dm = dm;
1303ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1304ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1305ccaff030SJeremy L Thompson   user->ceed = ceed;
1306ccaff030SJeremy L Thompson 
1307ea6e0f84SLeila Ghaffari   // Calculate qdataVol and ICs
1308ccaff030SJeremy L Thompson   // Set up state global and local vectors
1309ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1310ccaff030SJeremy L Thompson 
1311c96c872fSLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceedVol, Qloc); CHKERRQ(ierr);
1312ccaff030SJeremy L Thompson 
1313ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1314ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
1315ea6e0f84SLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdataVol, CEED_REQUEST_IMMEDIATE);
1316ea6e0f84SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictqVol, basisqVol, restrictqdiVol, qdataVol,
1317ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1318ccaff030SJeremy L Thompson 
1319c96c872fSLeila Ghaffari   ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceedVol, dm, Qloc, Q, restrictqVol,
13200c6c0b13SLeila Ghaffari                                &ctxSetup, 0.0);
1321ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1322ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1323ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1324ccaff030SJeremy L Thompson     // slower execution.
1325ccaff030SJeremy L Thompson     Vec Qbc;
1326ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1327ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1328ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1329ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1330ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1331ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1332ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
13330c6c0b13SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_NS); CHKERRQ(ierr);
1334ccaff030SJeremy L Thompson   }
1335ccaff030SJeremy L Thompson 
1336ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1337ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1338ccaff030SJeremy L Thompson   // Gather initial Q values
1339ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1340ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1341ccaff030SJeremy L Thompson     PetscViewer viewer;
1342ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1343ccaff030SJeremy L Thompson     // Read input
1344ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1345ccaff030SJeremy L Thompson                          user->outputfolder);
1346ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1347ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1348ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1349ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1350ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
13510c6c0b13SLeila Ghaffari   } else {
13520c6c0b13SLeila Ghaffari     //ierr = DMLocalToGlobal(dm, Qloc, INSERT_VALUES, Q);CHKERRQ(ierr);
1353ccaff030SJeremy L Thompson   }
1354ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1355ccaff030SJeremy L Thompson 
1356ccaff030SJeremy L Thompson   // Create and setup TS
1357ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1358ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1359ccaff030SJeremy L Thompson   if (implicit) {
1360ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1361ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1362ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1363ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1364ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1365ccaff030SJeremy L Thompson     }
1366ccaff030SJeremy L Thompson   } else {
1367ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1368ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1369ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1370ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1371ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1372ccaff030SJeremy L Thompson   }
1373ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1374ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1375ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
13760c6c0b13SLeila Ghaffari   if (test) {ierr = TSSetMaxSteps(ts, 1); CHKERRQ(ierr);}
1377ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1378ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1379ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1380ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1381ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1382ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
13830c6c0b13SLeila Ghaffari     if (!test) {
1384ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1385ccaff030SJeremy L Thompson     }
1386ccaff030SJeremy L Thompson   } else { // continue from time of last output
1387ccaff030SJeremy L Thompson     PetscReal time;
1388ccaff030SJeremy L Thompson     PetscInt count;
1389ccaff030SJeremy L Thompson     PetscViewer viewer;
1390ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1391ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1392ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1393ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1394ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1395ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1396ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1397ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1398ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1399ccaff030SJeremy L Thompson   }
14000c6c0b13SLeila Ghaffari   if (!test) {
1401ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1402ccaff030SJeremy L Thompson   }
1403ccaff030SJeremy L Thompson 
1404ccaff030SJeremy L Thompson   // Solve
1405ccaff030SJeremy L Thompson   start = MPI_Wtime();
1406ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1407ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1408ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1409ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts,&ftime); CHKERRQ(ierr);
1410ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1411ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
14120c6c0b13SLeila Ghaffari   if (!test) {
1413ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
14140c6c0b13SLeila Ghaffari                        "Time taken for solution: %g\n",
1415ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1416ccaff030SJeremy L Thompson   }
1417ccaff030SJeremy L Thompson 
1418ccaff030SJeremy L Thompson   // Get error
14190c6c0b13SLeila Ghaffari   if (problem->non_zero_time && !test) {
1420ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1421ccaff030SJeremy L Thompson     PetscReal norm;
1422ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1423ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1424ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1425ccaff030SJeremy L Thompson 
1426c96c872fSLeila Ghaffari     ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceedVol, dm, Qexactloc, Qexact,
1427ea6e0f84SLeila Ghaffari                                  restrictqVol, &ctxSetup, ftime); CHKERRQ(ierr);
1428ccaff030SJeremy L Thompson 
1429ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1430ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1431c96c872fSLeila Ghaffari     CeedVectorDestroy(&q0ceedVol);
1432ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1433ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1434ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
1435ccaff030SJeremy L Thompson   }
1436ccaff030SJeremy L Thompson 
1437ccaff030SJeremy L Thompson   // Output Statistics
1438ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
14390c6c0b13SLeila Ghaffari   if (!test) {
1440ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1441ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1442ccaff030SJeremy L Thompson                        steps,(double)ftime); CHKERRQ(ierr);
1443ccaff030SJeremy L Thompson   }
14449cf88b28Svaleriabarra 
1445ccaff030SJeremy L Thompson   // Clean up libCEED
1446ea6e0f84SLeila Ghaffari   CeedVectorDestroy(&qdataVol);
1447ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1448ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1449ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1450ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
1451ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisqVol);
1452ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxVol);
1453ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxcVol);
1454ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqVol);
1455ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictxVol);
1456ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdiVol);
1457ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1458ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1459ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1460ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1461ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1462ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
14636a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_vol);
14646a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_vol);
14656a0edaf9SLeila Ghaffari   CeedDestroy(&ceed);
14666a0edaf9SLeila Ghaffari 
14676a0edaf9SLeila Ghaffari   CeedVectorDestroy(&qdataSur);
14686a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisqSur);
14696a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxSur);
14706a0edaf9SLeila Ghaffari   CeedBasisDestroy(&basisxcSur);
14716a0edaf9SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqSur);
14726a0edaf9SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictxSur);
14736a0edaf9SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdiSur);
14746a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupSur);
14756a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsSur);
14766a0edaf9SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionSur);
14776a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&op_setupSur);
14786a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_rhs_sur);
14796a0edaf9SLeila Ghaffari   CeedOperatorDestroy(&user->op_ifunction_sur);
14806a0edaf9SLeila Ghaffari 
1481ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1482ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1483ccaff030SJeremy L Thompson 
1484ccaff030SJeremy L Thompson   // Clean up PETSc
1485ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1486ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1487ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1488ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1489ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1490ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1491ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1492ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1493ccaff030SJeremy L Thompson   return PetscFinalize();
1494ccaff030SJeremy L Thompson }
1495ccaff030SJeremy L Thompson 
1496