xref: /libCEED/examples/fluids/navierstokes.c (revision b0137797eb2ae643a21156d7b242c1fbfc1a0af9)
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"
47*b0137797SLeila 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 {
79ea6e0f84SLeila Ghaffari   CeedInt dim, qdatasizeVol, qdatasizeSur;
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,
93c96c872fSLeila Ghaffari     .qdatasizeSur              = 5,
94*b0137797SLeila Ghaffari     .setupVol                  = Setup,
95*b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
96*b0137797SLeila Ghaffari   //.setupSur                  = SetupBoundary,
97*b0137797SLeila Ghaffari   //.setupSur_loc              = SetupBoundary_loc,
98ccaff030SJeremy L Thompson     .ics                       = ICsDC,
99ccaff030SJeremy L Thompson     .ics_loc                   = ICsDC_loc,
100c96c872fSLeila Ghaffari     .applyVol_rhs              = DC,
101c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = DC_loc,
102ea6e0f84SLeila Ghaffari   //.applySur_rhs              = DC_Sur,
103ea6e0f84SLeila Ghaffari   //.applySur_rhs_loc          = DC_Sur_loc,
104c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_DC,
105c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_DC_loc,
106ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_DC_Sur,
107ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_DC_Sur_loc,
108ccaff030SJeremy L Thompson     .bc                        = Exact_DC,
1090c6c0b13SLeila Ghaffari     .non_zero_time             = false,
110ccaff030SJeremy L Thompson   },
111ccaff030SJeremy L Thompson   [NS_ADVECTION] = {
112ccaff030SJeremy L Thompson     .dim                       = 3,
113ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 10,
114c96c872fSLeila Ghaffari     .qdatasizeSur              = 5,
115*b0137797SLeila Ghaffari     .setupVol                  = Setup,
116*b0137797SLeila Ghaffari     .setupVol_loc              = Setup_loc,
117*b0137797SLeila Ghaffari   //.setupSur                  = SetupBoundary,
118*b0137797SLeila Ghaffari   //.setupSur_loc              = SetupBoundary_loc,
119ccaff030SJeremy L Thompson     .ics                       = ICsAdvection,
120ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection_loc,
121c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection,
122c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection_loc,
123ea6e0f84SLeila Ghaffari   //.applySur_rhs              = Advection_Sur,
124ea6e0f84SLeila Ghaffari   //.applySur_rhs_loc          = Advection_Sur_loc,
125c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection,
126c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection_loc,
127ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection_Sur,
128ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection_Sur_loc,
129ccaff030SJeremy L Thompson     .bc                        = Exact_Advection,
1300c6c0b13SLeila Ghaffari     .non_zero_time             = true,
131ccaff030SJeremy L Thompson   },
132ccaff030SJeremy L Thompson   [NS_ADVECTION2D] = {
133ccaff030SJeremy L Thompson     .dim                       = 2,
134ea6e0f84SLeila Ghaffari     .qdatasizeVol              = 5,
135*b0137797SLeila Ghaffari     .qdatasizeSur              = 1,
136c96c872fSLeila Ghaffari     .setupVol                  = Setup2d,
137c96c872fSLeila Ghaffari     .setupVol_loc              = Setup2d_loc,
138*b0137797SLeila Ghaffari     .setupSur                  = SetupBoundary2d,
139*b0137797SLeila Ghaffari     .setupSur_loc              = SetupBoundary2d_loc,
140ccaff030SJeremy L Thompson     .ics                       = ICsAdvection2d,
141ccaff030SJeremy L Thompson     .ics_loc                   = ICsAdvection2d_loc,
142c96c872fSLeila Ghaffari     .applyVol_rhs              = Advection2d,
143c96c872fSLeila Ghaffari     .applyVol_rhs_loc          = Advection2d_loc,
144*b0137797SLeila Ghaffari     .applySur_rhs              = Advection2d_Sur,
145*b0137797SLeila Ghaffari     .applySur_rhs_loc          = Advection2d_Sur_loc,
146c96c872fSLeila Ghaffari     .applyVol_ifunction        = IFunction_Advection2d,
147c96c872fSLeila Ghaffari     .applyVol_ifunction_loc    = IFunction_Advection2d_loc,
148ea6e0f84SLeila Ghaffari   //.applySur_ifunction        = IFunction_Advection2d_Sur,
149ea6e0f84SLeila Ghaffari   //.applySur_ifunction_loc    = IFunction_Advection2d_Sur_loc,
150ccaff030SJeremy L Thompson     .bc                        = Exact_Advection2d,
1510c6c0b13SLeila Ghaffari     .non_zero_time             = true,
152ccaff030SJeremy L Thompson   },
153ccaff030SJeremy L Thompson };
154ccaff030SJeremy L Thompson 
155ccaff030SJeremy L Thompson // PETSc user data
156ccaff030SJeremy L Thompson typedef struct User_ *User;
157ccaff030SJeremy L Thompson typedef struct Units_ *Units;
158ccaff030SJeremy L Thompson 
159ccaff030SJeremy L Thompson struct User_ {
160ccaff030SJeremy L Thompson   MPI_Comm comm;
161ccaff030SJeremy L Thompson   PetscInt outputfreq;
162ccaff030SJeremy L Thompson   DM dm;
163ccaff030SJeremy L Thompson   DM dmviz;
164ccaff030SJeremy L Thompson   Mat interpviz;
165ccaff030SJeremy L Thompson   Ceed ceed;
166ccaff030SJeremy L Thompson   Units units;
167ccaff030SJeremy L Thompson   CeedVector qceed, qdotceed, gceed;
168ea6e0f84SLeila Ghaffari   CeedOperator op_rhs_vol, op_rhs_sur, op_rhs,
169ea6e0f84SLeila Ghaffari                op_ifunction_vol, op_ifunction_sur, op_ifunction;
170ccaff030SJeremy L Thompson   Vec M;
171ccaff030SJeremy L Thompson   char outputfolder[PETSC_MAX_PATH_LEN];
172ccaff030SJeremy L Thompson   PetscInt contsteps;
173ccaff030SJeremy L Thompson };
174ccaff030SJeremy L Thompson 
175ccaff030SJeremy L Thompson struct Units_ {
176ccaff030SJeremy L Thompson   // fundamental units
177ccaff030SJeremy L Thompson   PetscScalar meter;
178ccaff030SJeremy L Thompson   PetscScalar kilogram;
179ccaff030SJeremy L Thompson   PetscScalar second;
180ccaff030SJeremy L Thompson   PetscScalar Kelvin;
181ccaff030SJeremy L Thompson   // derived units
182ccaff030SJeremy L Thompson   PetscScalar Pascal;
183ccaff030SJeremy L Thompson   PetscScalar JperkgK;
184ccaff030SJeremy L Thompson   PetscScalar mpersquareds;
185ccaff030SJeremy L Thompson   PetscScalar WpermK;
186ccaff030SJeremy L Thompson   PetscScalar kgpercubicm;
187ccaff030SJeremy L Thompson   PetscScalar kgpersquaredms;
188ccaff030SJeremy L Thompson   PetscScalar Joulepercubicm;
189ccaff030SJeremy L Thompson };
190ccaff030SJeremy L Thompson 
191ccaff030SJeremy L Thompson typedef struct SimpleBC_ *SimpleBC;
192ccaff030SJeremy L Thompson struct SimpleBC_ {
193ccaff030SJeremy L Thompson   PetscInt nwall, nslip[3];
1940c6c0b13SLeila Ghaffari   PetscInt walls[10], slips[3][10];
195ccaff030SJeremy L Thompson };
196ccaff030SJeremy L Thompson 
197ccaff030SJeremy L Thompson // Essential BC dofs are encoded in closure indices as -(i+1).
198ccaff030SJeremy L Thompson static PetscInt Involute(PetscInt i) {
199ccaff030SJeremy L Thompson   return i >= 0 ? i : -(i+1);
200ccaff030SJeremy L Thompson }
201ccaff030SJeremy L Thompson 
202ccaff030SJeremy L Thompson // Utility function to create local CEED restriction
203ccaff030SJeremy L Thompson static PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt P,
2040c6c0b13SLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value,
205ccaff030SJeremy L Thompson     CeedElemRestriction *Erestrict) {
206ccaff030SJeremy L Thompson 
207ccaff030SJeremy L Thompson   PetscSection   section;
2080c6c0b13SLeila Ghaffari   PetscInt       p, Nelem, Ndof, *erestrict, eoffset, nfields, dim,
2090c6c0b13SLeila Ghaffari                  depth;
2100c6c0b13SLeila Ghaffari   DMLabel depthLabel;
2110c6c0b13SLeila Ghaffari   IS depthIS, iterIS;
2120c6c0b13SLeila Ghaffari   const PetscInt *iterIndices;
213ccaff030SJeremy L Thompson   PetscErrorCode ierr;
214ccaff030SJeremy L Thompson   Vec Uloc;
215ccaff030SJeremy L Thompson 
216ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
217ccaff030SJeremy L Thompson   ierr = DMGetDimension(dm, &dim); CHKERRQ(ierr);
218ccaff030SJeremy L Thompson   ierr = DMGetLocalSection(dm,&section); CHKERRQ(ierr);
219ccaff030SJeremy L Thompson   ierr = PetscSectionGetNumFields(section, &nfields); CHKERRQ(ierr);
220ccaff030SJeremy L Thompson   PetscInt ncomp[nfields], fieldoff[nfields+1];
221ccaff030SJeremy L Thompson   fieldoff[0] = 0;
222ccaff030SJeremy L Thompson   for (PetscInt f=0; f<nfields; f++) {
223ccaff030SJeremy L Thompson     ierr = PetscSectionGetFieldComponents(section, f, &ncomp[f]); CHKERRQ(ierr);
224ccaff030SJeremy L Thompson     fieldoff[f+1] = fieldoff[f] + ncomp[f];
225ccaff030SJeremy L Thompson   }
226ccaff030SJeremy L Thompson 
2270c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepth(dm, &depth); CHKERRQ(ierr);
2280c6c0b13SLeila Ghaffari   ierr = DMPlexGetDepthLabel(dm, &depthLabel); CHKERRQ(ierr);
2290c6c0b13SLeila Ghaffari   ierr = DMLabelGetStratumIS(depthLabel, depth - height, &depthIS); CHKERRQ(ierr);
2300c6c0b13SLeila Ghaffari   if (domainLabel) {
2310c6c0b13SLeila Ghaffari     IS domainIS;
2320c6c0b13SLeila Ghaffari     ierr = DMLabelGetStratumIS(domainLabel, value, &domainIS); CHKERRQ(ierr);
2330c6c0b13SLeila Ghaffari     ierr = ISIntersect(depthIS, domainIS, &iterIS); CHKERRQ(ierr);
2340c6c0b13SLeila Ghaffari     ierr = ISDestroy(&domainIS); CHKERRQ(ierr);
2350c6c0b13SLeila Ghaffari     ierr = ISDestroy(&depthIS); CHKERRQ(ierr);
2360c6c0b13SLeila Ghaffari   } else {
2370c6c0b13SLeila Ghaffari     iterIS = depthIS;
2380c6c0b13SLeila Ghaffari   }
2390c6c0b13SLeila Ghaffari   ierr = ISGetLocalSize(iterIS, &Nelem); CHKERRQ(ierr);
2400c6c0b13SLeila Ghaffari   ierr = ISGetIndices(iterIS, &iterIndices); CHKERRQ(ierr);
241ccaff030SJeremy L Thompson   ierr = PetscMalloc1(Nelem*PetscPowInt(P, dim), &erestrict); CHKERRQ(ierr);
2420c6c0b13SLeila Ghaffari   for (p=0,eoffset=0; p<Nelem; p++) {
2430c6c0b13SLeila Ghaffari     PetscInt c = iterIndices[p];
244ccaff030SJeremy L Thompson     PetscInt numindices, *indices, nnodes;
2450c6c0b13SLeila Ghaffari     ierr = DMPlexGetClosureIndices(dm, section, section, c, &numindices,
2460c6c0b13SLeila Ghaffari                                    &indices, NULL); CHKERRQ(ierr);
2470c6c0b13SLeila Ghaffari     if (numindices % fieldoff[nfields])
2480c6c0b13SLeila Ghaffari       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP,
2490c6c0b13SLeila Ghaffari                "Number of closure indices not compatible with Cell %D", c);
250ccaff030SJeremy L Thompson     nnodes = numindices / fieldoff[nfields];
251ccaff030SJeremy L Thompson     for (PetscInt i=0; i<nnodes; i++) {
252ccaff030SJeremy L Thompson       // Check that indices are blocked by node and thus can be coalesced as a single field with
253ccaff030SJeremy L Thompson       // fieldoff[nfields] = sum(ncomp) components.
254ccaff030SJeremy L Thompson       for (PetscInt f=0; f<nfields; f++) {
255ccaff030SJeremy L Thompson         for (PetscInt j=0; j<ncomp[f]; j++) {
256ccaff030SJeremy L Thompson           if (Involute(indices[fieldoff[f]*nnodes + i*ncomp[f] + j])
257ccaff030SJeremy L Thompson               != Involute(indices[i*ncomp[0]]) + fieldoff[f] + j)
258ccaff030SJeremy L Thompson             SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
259ccaff030SJeremy L Thompson                      "Cell %D closure indices not interlaced for node %D field %D component %D",
260ccaff030SJeremy L Thompson                      c, i, f, j);
261ccaff030SJeremy L Thompson         }
262ccaff030SJeremy L Thompson       }
263ccaff030SJeremy L Thompson       // Essential boundary conditions are encoded as -(loc+1), but we don't care so we decode.
264ccaff030SJeremy L Thompson       PetscInt loc = Involute(indices[i*ncomp[0]]);
2656f55dfd5Svaleriabarra       erestrict[eoffset++] = loc;
266ccaff030SJeremy L Thompson     }
2670c6c0b13SLeila Ghaffari     ierr = DMPlexRestoreClosureIndices(dm, section, section, c, &numindices,
2680c6c0b13SLeila Ghaffari                                        &indices, NULL); CHKERRQ(ierr);
269ccaff030SJeremy L Thompson   }
2700c6c0b13SLeila Ghaffari   if (eoffset != Nelem*PetscPowInt(P, dim))
2710c6c0b13SLeila Ghaffari     SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_LIB,
2720c6c0b13SLeila Ghaffari              "ElemRestriction of size (%D,%D) initialized %D nodes", Nelem,
273ccaff030SJeremy L Thompson              PetscPowInt(P, dim),eoffset);
2740c6c0b13SLeila Ghaffari   ierr = ISRestoreIndices(iterIS, &iterIndices); CHKERRQ(ierr);
2750c6c0b13SLeila Ghaffari   ierr = ISDestroy(&iterIS); CHKERRQ(ierr);
2760c6c0b13SLeila Ghaffari 
277ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Uloc); CHKERRQ(ierr);
278ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(Uloc, &Ndof); CHKERRQ(ierr);
279ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Uloc); CHKERRQ(ierr);
2806f55dfd5Svaleriabarra   CeedElemRestrictionCreate(ceed, Nelem, PetscPowInt(P, dim), fieldoff[nfields],
2816f55dfd5Svaleriabarra                             1, Ndof, CEED_MEM_HOST, CEED_COPY_VALUES, erestrict,
2826f55dfd5Svaleriabarra                             Erestrict);
283ccaff030SJeremy L Thompson   ierr = PetscFree(erestrict); CHKERRQ(ierr);
284ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
285ccaff030SJeremy L Thompson }
286ccaff030SJeremy L Thompson 
287c96c872fSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
288c96c872fSLeila Ghaffari static PetscErrorCode GetRestriction(Ceed ceed, DM dm, CeedInt ncompx, CeedInt dim,
289c96c872fSLeila Ghaffari     CeedInt height, DMLabel domainLabel, CeedInt value, CeedInt P, CeedInt Q,
290c96c872fSLeila Ghaffari     CeedInt qdatasize, CeedElemRestriction *restrictq,
291c96c872fSLeila Ghaffari     CeedElemRestriction *restrictx, CeedElemRestriction *restrictqdi) {
292c96c872fSLeila Ghaffari 
293c96c872fSLeila Ghaffari   DM dmcoord;
294c96c872fSLeila Ghaffari   CeedInt localNelem;
295c96c872fSLeila Ghaffari   CeedElemRestriction  *restrictxcoord = NULL;
296c96c872fSLeila Ghaffari   CeedInt Qdim = CeedIntPow(Q, dim);
297c96c872fSLeila Ghaffari   //CeedInt numPdim = CeedIntPow(P, dim);
298c96c872fSLeila Ghaffari   //Vec Xloc = NULL;
299c96c872fSLeila Ghaffari   PetscErrorCode ierr;
300c96c872fSLeila Ghaffari 
301c96c872fSLeila Ghaffari   PetscFunctionBeginUser;
302c96c872fSLeila Ghaffari   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
303c96c872fSLeila Ghaffari   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
304c96c872fSLeila Ghaffari   CHKERRQ(ierr);
305c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dm, P, height, domainLabel, value, restrictq);
306c96c872fSLeila Ghaffari   CHKERRQ(ierr);
307c96c872fSLeila Ghaffari   ierr = CreateRestrictionFromPlex(ceed, dmcoord, 2, height, domainLabel, value, restrictx);
308c96c872fSLeila Ghaffari   CHKERRQ(ierr);
309c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(*restrictq, &localNelem);
310c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, Qdim,
311c96c872fSLeila Ghaffari                                    qdatasize, qdatasize*localNelem*Qdim,
312c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictqdi);
313c96c872fSLeila Ghaffari   CeedElemRestrictionCreateStrided(ceed, localNelem, PetscPowInt(P, dim),
314c96c872fSLeila Ghaffari                                    ncompx,
315c96c872fSLeila Ghaffari                                    ncompx*localNelem*PetscPowInt(P, dim),
316c96c872fSLeila Ghaffari                                    CEED_STRIDES_BACKEND, restrictxcoord);
317c96c872fSLeila Ghaffari   PetscFunctionReturn(0);
318c96c872fSLeila Ghaffari }
319c96c872fSLeila Ghaffari 
320ccaff030SJeremy L Thompson static int CreateVectorFromPetscVec(Ceed ceed, Vec p, CeedVector *v) {
321ccaff030SJeremy L Thompson   PetscErrorCode ierr;
322ccaff030SJeremy L Thompson   PetscInt m;
323ccaff030SJeremy L Thompson 
324ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
325ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &m); CHKERRQ(ierr);
326ccaff030SJeremy L Thompson   ierr = CeedVectorCreate(ceed, m, v); CHKERRQ(ierr);
327ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
328ccaff030SJeremy L Thompson }
329ccaff030SJeremy L Thompson 
330ccaff030SJeremy L Thompson static int VectorPlacePetscVec(CeedVector c, Vec p) {
331ccaff030SJeremy L Thompson   PetscErrorCode ierr;
332ccaff030SJeremy L Thompson   PetscInt mceed,mpetsc;
333ccaff030SJeremy L Thompson   PetscScalar *a;
334ccaff030SJeremy L Thompson 
335ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
336ccaff030SJeremy L Thompson   ierr = CeedVectorGetLength(c, &mceed); CHKERRQ(ierr);
337ccaff030SJeremy L Thompson   ierr = VecGetLocalSize(p, &mpetsc); CHKERRQ(ierr);
338ccaff030SJeremy L Thompson   if (mceed != mpetsc) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,
3390c6c0b13SLeila Ghaffari                                   "Cannot place PETSc Vec of length %D in CeedVector of length %D",mpetsc,mceed);
340ccaff030SJeremy L Thompson   ierr = VecGetArray(p, &a); CHKERRQ(ierr);
341ccaff030SJeremy L Thompson   CeedVectorSetArray(c, CEED_MEM_HOST, CEED_USE_POINTER, a);
342ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
343ccaff030SJeremy L Thompson }
344ccaff030SJeremy L Thompson 
345ccaff030SJeremy L Thompson static PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm,
346ccaff030SJeremy L Thompson     PetscBool insertEssential, Vec Qloc, PetscReal time, Vec faceGeomFVM,
347ccaff030SJeremy L Thompson     Vec cellGeomFVM, Vec gradFVM) {
348ccaff030SJeremy L Thompson   PetscErrorCode ierr;
349ccaff030SJeremy L Thompson   Vec Qbc;
350ccaff030SJeremy L Thompson 
351ccaff030SJeremy L Thompson   PetscFunctionBegin;
352ccaff030SJeremy L Thompson   ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
353ccaff030SJeremy L Thompson   ierr = VecAXPY(Qloc, 1., Qbc); CHKERRQ(ierr);
354ccaff030SJeremy L Thompson   ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
355ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
356ccaff030SJeremy L Thompson }
357ccaff030SJeremy L Thompson 
358ccaff030SJeremy L Thompson // This is the RHS of the ODE, given as u_t = G(t,u)
359ccaff030SJeremy L Thompson // This function takes in a state vector Q and writes into G
360ccaff030SJeremy L Thompson static PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *userData) {
361ccaff030SJeremy L Thompson   PetscErrorCode ierr;
362ccaff030SJeremy L Thompson   User user = *(User *)userData;
363ccaff030SJeremy L Thompson   PetscScalar *q, *g;
364ccaff030SJeremy L Thompson   Vec Qloc, Gloc;
365ccaff030SJeremy L Thompson 
366ccaff030SJeremy L Thompson   // Global-to-local
367ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
368ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
369ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
370ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
371ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
372ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
373ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
374ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
375ccaff030SJeremy L Thompson 
376ccaff030SJeremy L Thompson   // Ceed Vectors
377ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
378ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
379ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER, q);
380ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
381ccaff030SJeremy L Thompson 
382ccaff030SJeremy L Thompson   // Apply CEED operator
383ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_rhs, user->qceed, user->gceed,
384ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
385ccaff030SJeremy L Thompson 
386ccaff030SJeremy L Thompson   // Restore vectors
387ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, (const PetscScalar **)&q); CHKERRQ(ierr);
388ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
389ccaff030SJeremy L Thompson 
390ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
391ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
392ccaff030SJeremy L Thompson 
393ccaff030SJeremy L Thompson   // Inverse of the lumped mass matrix
394ccaff030SJeremy L Thompson   ierr = VecPointwiseMult(G, G, user->M); // M is Minv
395ccaff030SJeremy L Thompson   CHKERRQ(ierr);
396ccaff030SJeremy L Thompson 
397ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
398ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
399ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
400ccaff030SJeremy L Thompson }
401ccaff030SJeremy L Thompson 
402ccaff030SJeremy L Thompson static PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Qdot, Vec G,
403ccaff030SJeremy L Thompson                                    void *userData) {
404ccaff030SJeremy L Thompson   PetscErrorCode ierr;
405ccaff030SJeremy L Thompson   User user = *(User *)userData;
406ccaff030SJeremy L Thompson   const PetscScalar *q, *qdot;
407ccaff030SJeremy L Thompson   PetscScalar *g;
408ccaff030SJeremy L Thompson   Vec Qloc, Qdotloc, Gloc;
409ccaff030SJeremy L Thompson 
410ccaff030SJeremy L Thompson   // Global-to-local
411ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
412ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
413ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
414ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
415ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
416ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
417ccaff030SJeremy L Thompson   ierr = DMPlexInsertBoundaryValues(user->dm, PETSC_TRUE, Qloc, 0.0,
418ccaff030SJeremy L Thompson                                     NULL, NULL, NULL); CHKERRQ(ierr);
419ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qdotloc); CHKERRQ(ierr);
420ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Qdot, INSERT_VALUES, Qdotloc); CHKERRQ(ierr);
421ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Gloc); CHKERRQ(ierr);
422ccaff030SJeremy L Thompson 
423ccaff030SJeremy L Thompson   // Ceed Vectors
424ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qloc, &q); CHKERRQ(ierr);
425ccaff030SJeremy L Thompson   ierr = VecGetArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
426ccaff030SJeremy L Thompson   ierr = VecGetArray(Gloc, &g); CHKERRQ(ierr);
427ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qceed, CEED_MEM_HOST, CEED_USE_POINTER,
428ccaff030SJeremy L Thompson                      (PetscScalar *)q);
429ccaff030SJeremy L Thompson   CeedVectorSetArray(user->qdotceed, CEED_MEM_HOST, CEED_USE_POINTER,
430ccaff030SJeremy L Thompson                      (PetscScalar *)qdot);
431ccaff030SJeremy L Thompson   CeedVectorSetArray(user->gceed, CEED_MEM_HOST, CEED_USE_POINTER, g);
432ccaff030SJeremy L Thompson 
433ccaff030SJeremy L Thompson   // Apply CEED operator
434ccaff030SJeremy L Thompson   CeedOperatorApply(user->op_ifunction, user->qceed, user->gceed,
435ccaff030SJeremy L Thompson                     CEED_REQUEST_IMMEDIATE);
436ccaff030SJeremy L Thompson 
437ccaff030SJeremy L Thompson   // Restore vectors
438ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qloc, &q); CHKERRQ(ierr);
439ccaff030SJeremy L Thompson   ierr = VecRestoreArrayRead(Qdotloc, &qdot); CHKERRQ(ierr);
440ccaff030SJeremy L Thompson   ierr = VecRestoreArray(Gloc, &g); CHKERRQ(ierr);
441ccaff030SJeremy L Thompson 
442ccaff030SJeremy L Thompson   ierr = VecZeroEntries(G); CHKERRQ(ierr);
443ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(user->dm, Gloc, ADD_VALUES, G); CHKERRQ(ierr);
444ccaff030SJeremy L Thompson 
445ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
446ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qdotloc); CHKERRQ(ierr);
447ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Gloc); CHKERRQ(ierr);
448ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
449ccaff030SJeremy L Thompson }
450ccaff030SJeremy L Thompson 
451ccaff030SJeremy L Thompson // User provided TS Monitor
452ccaff030SJeremy L Thompson static PetscErrorCode TSMonitor_NS(TS ts, PetscInt stepno, PetscReal time,
453ccaff030SJeremy L Thompson                                    Vec Q, void *ctx) {
454ccaff030SJeremy L Thompson   User user = ctx;
455ccaff030SJeremy L Thompson   Vec Qloc;
456ccaff030SJeremy L Thompson   char filepath[PETSC_MAX_PATH_LEN];
457ccaff030SJeremy L Thompson   PetscViewer viewer;
458ccaff030SJeremy L Thompson   PetscErrorCode ierr;
459ccaff030SJeremy L Thompson 
460ccaff030SJeremy L Thompson   // Set up output
461ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
462ccaff030SJeremy L Thompson   // Print every 'outputfreq' steps
463ccaff030SJeremy L Thompson   if (stepno % user->outputfreq != 0)
464ccaff030SJeremy L Thompson     PetscFunctionReturn(0);
465ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
466ccaff030SJeremy L Thompson   ierr = PetscObjectSetName((PetscObject)Qloc, "StateVec"); CHKERRQ(ierr);
467ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
468ccaff030SJeremy L Thompson   ierr = DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
469ccaff030SJeremy L Thompson 
470ccaff030SJeremy L Thompson   // Output
471ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-%03D.vtu",
472ccaff030SJeremy L Thompson                        user->outputfolder, stepno + user->contsteps);
473ccaff030SJeremy L Thompson   CHKERRQ(ierr);
474ccaff030SJeremy L Thompson   ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Q), filepath,
475ccaff030SJeremy L Thompson                             FILE_MODE_WRITE, &viewer); CHKERRQ(ierr);
476ccaff030SJeremy L Thompson   ierr = VecView(Qloc, viewer); CHKERRQ(ierr);
477ccaff030SJeremy L Thompson   if (user->dmviz) {
478ccaff030SJeremy L Thompson     Vec Qrefined, Qrefined_loc;
479ccaff030SJeremy L Thompson     char filepath_refined[PETSC_MAX_PATH_LEN];
480ccaff030SJeremy L Thompson     PetscViewer viewer_refined;
481ccaff030SJeremy L Thompson 
482ccaff030SJeremy L Thompson     ierr = DMGetGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
483ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
484ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)Qrefined_loc, "Refined");
485ccaff030SJeremy L Thompson     CHKERRQ(ierr);
486ccaff030SJeremy L Thompson     ierr = MatInterpolate(user->interpviz, Q, Qrefined); CHKERRQ(ierr);
487ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qrefined_loc); CHKERRQ(ierr);
488ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(user->dmviz, Qrefined, INSERT_VALUES, Qrefined_loc);
489ccaff030SJeremy L Thompson     CHKERRQ(ierr);
490ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath_refined, sizeof filepath_refined,
491ccaff030SJeremy L Thompson                          "%s/nsrefined-%03D.vtu",
492ccaff030SJeremy L Thompson                          user->outputfolder, stepno + user->contsteps);
493ccaff030SJeremy L Thompson     CHKERRQ(ierr);
494ccaff030SJeremy L Thompson     ierr = PetscViewerVTKOpen(PetscObjectComm((PetscObject)Qrefined),
495ccaff030SJeremy L Thompson                               filepath_refined,
496ccaff030SJeremy L Thompson                               FILE_MODE_WRITE, &viewer_refined); CHKERRQ(ierr);
497ccaff030SJeremy L Thompson     ierr = VecView(Qrefined_loc, viewer_refined); CHKERRQ(ierr);
498ccaff030SJeremy L Thompson     ierr = DMRestoreLocalVector(user->dmviz, &Qrefined_loc); CHKERRQ(ierr);
499ccaff030SJeremy L Thompson     ierr = DMRestoreGlobalVector(user->dmviz, &Qrefined); CHKERRQ(ierr);
500ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer_refined); CHKERRQ(ierr);
501ccaff030SJeremy L Thompson   }
5020c6c0b13SLeila Ghaffari   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
503ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(user->dm, &Qloc); CHKERRQ(ierr);
504ccaff030SJeremy L Thompson 
505ccaff030SJeremy L Thompson   // Save data in a binary file for continuation of simulations
506ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.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   ierr = VecView(Q, viewer); CHKERRQ(ierr);
511ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
512ccaff030SJeremy L Thompson 
513ccaff030SJeremy L Thompson   // Save time stamp
514ccaff030SJeremy L Thompson   // Dimensionalize time back
515ccaff030SJeremy L Thompson   time /= user->units->second;
516ccaff030SJeremy L Thompson   ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
517ccaff030SJeremy L Thompson                        user->outputfolder); CHKERRQ(ierr);
518ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryOpen(user->comm, filepath, FILE_MODE_WRITE, &viewer);
519ccaff030SJeremy L Thompson   CHKERRQ(ierr);
520ccaff030SJeremy L Thompson   #if PETSC_VERSION_GE(3,13,0)
521ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL);
522ccaff030SJeremy L Thompson   #else
523ccaff030SJeremy L Thompson   ierr = PetscViewerBinaryWrite(viewer, &time, 1, PETSC_REAL, true);
524ccaff030SJeremy L Thompson   #endif
525ccaff030SJeremy L Thompson   CHKERRQ(ierr);
526ccaff030SJeremy L Thompson   ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
527ccaff030SJeremy L Thompson 
528ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
529ccaff030SJeremy L Thompson }
530ccaff030SJeremy L Thompson 
5310c6c0b13SLeila Ghaffari static PetscErrorCode ICs_PetscMultiplicity(CeedOperator op_ics,
532ccaff030SJeremy L Thompson     CeedVector xcorners, CeedVector q0ceed, DM dm, Vec Qloc, Vec Q,
533ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, SetupContext ctxSetup, CeedScalar time) {
534ccaff030SJeremy L Thompson   PetscErrorCode ierr;
535ccaff030SJeremy L Thompson   CeedVector multlvec;
536ccaff030SJeremy L Thompson   Vec Multiplicity, MultiplicityLoc;
537ccaff030SJeremy L Thompson 
538ccaff030SJeremy L Thompson   ctxSetup->time = time;
539ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
540ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(q0ceed, Qloc); CHKERRQ(ierr);
541ccaff030SJeremy L Thompson   CeedOperatorApply(op_ics, xcorners, q0ceed, CEED_REQUEST_IMMEDIATE);
542ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
543ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Qloc, ADD_VALUES, Q); CHKERRQ(ierr);
544ccaff030SJeremy L Thompson 
545ccaff030SJeremy L Thompson   // Fix multiplicity for output of ICs
546ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
547ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &multlvec, NULL);
548ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(multlvec, MultiplicityLoc); CHKERRQ(ierr);
549ccaff030SJeremy L Thompson   CeedElemRestrictionGetMultiplicity(restrictq, multlvec);
550ccaff030SJeremy L Thompson   CeedVectorDestroy(&multlvec);
551ccaff030SJeremy L Thompson   ierr = DMGetGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
552ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Multiplicity); CHKERRQ(ierr);
553ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, MultiplicityLoc, ADD_VALUES, Multiplicity);
554ccaff030SJeremy L Thompson   CHKERRQ(ierr);
555ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Q, Q, Multiplicity); CHKERRQ(ierr);
556ccaff030SJeremy L Thompson   ierr = VecPointwiseDivide(Qloc, Qloc, MultiplicityLoc); CHKERRQ(ierr);
557ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &MultiplicityLoc); CHKERRQ(ierr);
558ccaff030SJeremy L Thompson   ierr = DMRestoreGlobalVector(dm, &Multiplicity); CHKERRQ(ierr);
559ccaff030SJeremy L Thompson 
560ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
561ccaff030SJeremy L Thompson }
562ccaff030SJeremy L Thompson 
563ccaff030SJeremy L Thompson static PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm,
564ccaff030SJeremy L Thompson     CeedElemRestriction restrictq, CeedBasis basisq,
565ccaff030SJeremy L Thompson     CeedElemRestriction restrictqdi, CeedVector qdata, Vec M) {
566ccaff030SJeremy L Thompson   PetscErrorCode ierr;
567ccaff030SJeremy L Thompson   CeedQFunction qf_mass;
568ccaff030SJeremy L Thompson   CeedOperator op_mass;
569ccaff030SJeremy L Thompson   CeedVector mceed;
570ccaff030SJeremy L Thompson   Vec Mloc;
571ccaff030SJeremy L Thompson   CeedInt ncompq, qdatasize;
572ccaff030SJeremy L Thompson 
573ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
574ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictq, &ncompq);
575ccaff030SJeremy L Thompson   CeedElemRestrictionGetNumComponents(restrictqdi, &qdatasize);
576ccaff030SJeremy L Thompson   // Create the Q-function that defines the action of the mass operator
577ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, Mass, Mass_loc, &qf_mass);
578ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "q", ncompq, CEED_EVAL_INTERP);
579ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_mass, "qdata", qdatasize, CEED_EVAL_NONE);
580ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_mass, "v", ncompq, CEED_EVAL_INTERP);
581ccaff030SJeremy L Thompson 
582ccaff030SJeremy L Thompson   // Create the mass operator
583ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
584ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "q", restrictq, basisq, CEED_VECTOR_ACTIVE);
585ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "qdata", restrictqdi,
586ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, qdata);
587ccaff030SJeremy L Thompson   CeedOperatorSetField(op_mass, "v", restrictq, basisq, CEED_VECTOR_ACTIVE);
588ccaff030SJeremy L Thompson 
589ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Mloc); CHKERRQ(ierr);
590ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Mloc); CHKERRQ(ierr);
591ccaff030SJeremy L Thompson   CeedElemRestrictionCreateVector(restrictq, &mceed, NULL);
592ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(mceed, Mloc); CHKERRQ(ierr);
593ccaff030SJeremy L Thompson 
594ccaff030SJeremy L Thompson   {
595ccaff030SJeremy L Thompson     // Compute a lumped mass matrix
596ccaff030SJeremy L Thompson     CeedVector onesvec;
597ccaff030SJeremy L Thompson     CeedElemRestrictionCreateVector(restrictq, &onesvec, NULL);
598ccaff030SJeremy L Thompson     CeedVectorSetValue(onesvec, 1.0);
599ccaff030SJeremy L Thompson     CeedOperatorApply(op_mass, onesvec, mceed, CEED_REQUEST_IMMEDIATE);
600ccaff030SJeremy L Thompson     CeedVectorDestroy(&onesvec);
601ccaff030SJeremy L Thompson     CeedOperatorDestroy(&op_mass);
602ccaff030SJeremy L Thompson     CeedVectorDestroy(&mceed);
603ccaff030SJeremy L Thompson   }
604ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_mass);
605ccaff030SJeremy L Thompson 
606ccaff030SJeremy L Thompson   ierr = VecZeroEntries(M); CHKERRQ(ierr);
607ccaff030SJeremy L Thompson   ierr = DMLocalToGlobal(dm, Mloc, ADD_VALUES, M); CHKERRQ(ierr);
608ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Mloc); CHKERRQ(ierr);
609ccaff030SJeremy L Thompson 
610ccaff030SJeremy L Thompson   // Invert diagonally lumped mass vector for RHS function
611ccaff030SJeremy L Thompson   ierr = VecReciprocal(M); CHKERRQ(ierr);
612ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
613ccaff030SJeremy L Thompson }
614ccaff030SJeremy L Thompson 
6150c6c0b13SLeila Ghaffari PetscErrorCode SetUpDM(DM dm, problemData *problem, PetscInt degree,
616ff6701fcSJed Brown                        SimpleBC bc, void *ctxSetup) {
617ccaff030SJeremy L Thompson   PetscErrorCode ierr;
618ccaff030SJeremy L Thompson 
619ccaff030SJeremy L Thompson   PetscFunctionBeginUser;
620ccaff030SJeremy L Thompson   {
621ccaff030SJeremy L Thompson     // Configure the finite element space and boundary conditions
622ccaff030SJeremy L Thompson     PetscFE fe;
6230c6c0b13SLeila Ghaffari     PetscSpace fespace;
624ccaff030SJeremy L Thompson     PetscInt ncompq = 5;
625ff6701fcSJed Brown     ierr = PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim, ncompq,
626ff6701fcSJed Brown                                  PETSC_FALSE, degree, PETSC_DECIDE,
62732ed2d11SJed Brown                                  &fe); CHKERRQ(ierr);
628ccaff030SJeremy L Thompson     ierr = PetscObjectSetName((PetscObject)fe, "Q"); CHKERRQ(ierr);
629ccaff030SJeremy L Thompson     ierr = DMAddField(dm,NULL,(PetscObject)fe); CHKERRQ(ierr);
630ccaff030SJeremy L Thompson     ierr = DMCreateDS(dm); CHKERRQ(ierr);
6310c6c0b13SLeila Ghaffari     /* Wall boundary conditions are zero velocity and zero flux for density and energy */
6320c6c0b13SLeila Ghaffari     {
6330c6c0b13SLeila Ghaffari       PetscInt comps[3] = {1, 2, 3};
6340c6c0b13SLeila Ghaffari       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", "Face Sets", 0,
6350c6c0b13SLeila Ghaffari                            3, comps, (void(*)(void))problem->bc,
6360c6c0b13SLeila Ghaffari                            bc->nwall, bc->walls, ctxSetup); CHKERRQ(ierr);
6370c6c0b13SLeila Ghaffari     }
63807af6069Svaleriabarra     {
63907af6069Svaleriabarra       PetscInt comps[1] = {1};
64007af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipx", "Face Sets", 0,
64107af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[0],
64207af6069Svaleriabarra                            bc->slips[0], ctxSetup); CHKERRQ(ierr);
64307af6069Svaleriabarra       comps[0] = 2;
64407af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipy", "Face Sets", 0,
64507af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[1],
64607af6069Svaleriabarra                            bc->slips[1], ctxSetup); CHKERRQ(ierr);
64707af6069Svaleriabarra       comps[0] = 3;
64807af6069Svaleriabarra       ierr = DMAddBoundary(dm, DM_BC_ESSENTIAL, "slipz", "Face Sets", 0,
64907af6069Svaleriabarra                            1, comps, (void(*)(void))NULL, bc->nslip[2],
65007af6069Svaleriabarra                            bc->slips[2], ctxSetup); CHKERRQ(ierr);
65107af6069Svaleriabarra     }
652ccaff030SJeremy L Thompson     ierr = DMPlexSetClosurePermutationTensor(dm, PETSC_DETERMINE,NULL);
653ccaff030SJeremy L Thompson     CHKERRQ(ierr);
6540c6c0b13SLeila Ghaffari     ierr = PetscFEGetBasisSpace(fe, &fespace); CHKERRQ(ierr);
655ccaff030SJeremy L Thompson     ierr = PetscFEDestroy(&fe); CHKERRQ(ierr);
656ccaff030SJeremy L Thompson   }
657ccaff030SJeremy L Thompson   {
658ccaff030SJeremy L Thompson     // Empty name for conserved field (because there is only one field)
659ccaff030SJeremy L Thompson     PetscSection section;
660ccaff030SJeremy L Thompson     ierr = DMGetLocalSection(dm, &section); CHKERRQ(ierr);
661ccaff030SJeremy L Thompson     ierr = PetscSectionSetFieldName(section, 0, ""); CHKERRQ(ierr);
662ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 0, "Density");
663ccaff030SJeremy L Thompson     CHKERRQ(ierr);
664ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 1, "MomentumX");
665ccaff030SJeremy L Thompson     CHKERRQ(ierr);
666ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 2, "MomentumY");
667ccaff030SJeremy L Thompson     CHKERRQ(ierr);
668ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 3, "MomentumZ");
669ccaff030SJeremy L Thompson     CHKERRQ(ierr);
670ccaff030SJeremy L Thompson     ierr = PetscSectionSetComponentName(section, 0, 4, "EnergyDensity");
671ccaff030SJeremy L Thompson     CHKERRQ(ierr);
672ccaff030SJeremy L Thompson   }
673ccaff030SJeremy L Thompson   PetscFunctionReturn(0);
674ccaff030SJeremy L Thompson }
675ccaff030SJeremy L Thompson 
676ccaff030SJeremy L Thompson int main(int argc, char **argv) {
677ccaff030SJeremy L Thompson   PetscInt ierr;
678ccaff030SJeremy L Thompson   MPI_Comm comm;
6790c6c0b13SLeila Ghaffari   DM dm, dmcoord, dmviz, dmvizfine;
680ccaff030SJeremy L Thompson   Mat interpviz;
681ccaff030SJeremy L Thompson   TS ts;
682ccaff030SJeremy L Thompson   TSAdapt adapt;
683ccaff030SJeremy L Thompson   User user;
684ccaff030SJeremy L Thompson   Units units;
685ccaff030SJeremy L Thompson   char ceedresource[4096] = "/cpu/self";
686c96c872fSLeila Ghaffari   PetscInt localNelemVol, localNelemSur, lnodes, steps;
687ccaff030SJeremy L Thompson   const PetscInt ncompq = 5;
688ccaff030SJeremy L Thompson   PetscMPIInt rank;
689ccaff030SJeremy L Thompson   PetscScalar ftime;
690ccaff030SJeremy L Thompson   Vec Q, Qloc, Xloc;
691ccaff030SJeremy L Thompson   Ceed ceed;
692ea6e0f84SLeila Ghaffari   CeedInt numP_Vol, numP_Sur, numQ_Vol, numQ_Sur;
693c96c872fSLeila Ghaffari   CeedVector xcorners, qdataVol, qdataSur, q0ceedVol, q0ceedSur;
694ea6e0f84SLeila Ghaffari   CeedBasis basisxVol, basisxcVol, basisqVol, basisxSur, basisxcSur, basisqSur;
695ea6e0f84SLeila Ghaffari   CeedElemRestriction restrictxVol, restrictxcoordVol, restrictqVol, restrictqdiVol,
696ea6e0f84SLeila Ghaffari                       restrictxSur, restrictxcoordSur, restrictqSur, restrictqdiSur;
697ea6e0f84SLeila Ghaffari   CeedQFunction qf_setupVol, qf_setupSur, qf_ics, qf_rhsVol, qf_rhsSur,
698ea6e0f84SLeila Ghaffari                 qf_ifunctionVol, qf_ifunctionSur;
699ea6e0f84SLeila Ghaffari   CeedOperator op_setupVol, op_setupSur, op_ics;
700ccaff030SJeremy L Thompson   CeedScalar Rd;
701ccaff030SJeremy L Thompson   PetscScalar WpermK, Pascal, JperkgK, mpersquareds, kgpercubicm,
702ccaff030SJeremy L Thompson               kgpersquaredms, Joulepercubicm;
703ccaff030SJeremy L Thompson   problemType problemChoice;
704ccaff030SJeremy L Thompson   problemData *problem = NULL;
705ccaff030SJeremy L Thompson   StabilizationType stab;
7060c6c0b13SLeila Ghaffari   PetscBool   test, implicit;
707cb3e2689Svaleriabarra   PetscInt    viz_refine = 0;
708ccaff030SJeremy L Thompson   struct SimpleBC_ bc = {
7090c6c0b13SLeila Ghaffari     .nwall = 6,
7100c6c0b13SLeila Ghaffari     .walls = {1,2,3,4,5,6},
711ccaff030SJeremy L Thompson   };
712ccaff030SJeremy L Thompson   double start, cpu_time_used;
713ccaff030SJeremy L Thompson 
714ccaff030SJeremy L Thompson   // Create the libCEED contexts
715ccaff030SJeremy L Thompson   PetscScalar meter     = 1e-2;     // 1 meter in scaled length units
716ccaff030SJeremy L Thompson   PetscScalar second    = 1e-2;     // 1 second in scaled time units
717ccaff030SJeremy L Thompson   PetscScalar kilogram  = 1e-6;     // 1 kilogram in scaled mass units
718ccaff030SJeremy L Thompson   PetscScalar Kelvin    = 1;        // 1 Kelvin in scaled temperature units
719ccaff030SJeremy L Thompson   CeedScalar theta0     = 300.;     // K
720ccaff030SJeremy L Thompson   CeedScalar thetaC     = -15.;     // K
721ccaff030SJeremy L Thompson   CeedScalar P0         = 1.e5;     // Pa
722ccaff030SJeremy L Thompson   CeedScalar N          = 0.01;     // 1/s
723ccaff030SJeremy L Thompson   CeedScalar cv         = 717.;     // J/(kg K)
724ccaff030SJeremy L Thompson   CeedScalar cp         = 1004.;    // J/(kg K)
725ccaff030SJeremy L Thompson   CeedScalar g          = 9.81;     // m/s^2
726ccaff030SJeremy L Thompson   CeedScalar lambda     = -2./3.;   // -
727ccaff030SJeremy L Thompson   CeedScalar mu         = 75.;      // Pa s, dynamic viscosity
728ccaff030SJeremy L Thompson   // mu = 75 is not physical for air, but is good for numerical stability
729ccaff030SJeremy L Thompson   CeedScalar k          = 0.02638;  // W/(m K)
730ccaff030SJeremy L Thompson   CeedScalar CtauS      = 0.;       // dimensionless
731ccaff030SJeremy L Thompson   CeedScalar strong_form = 0.;      // [0,1]
732ccaff030SJeremy L Thompson   PetscScalar lx        = 8000.;    // m
733ccaff030SJeremy L Thompson   PetscScalar ly        = 8000.;    // m
734ccaff030SJeremy L Thompson   PetscScalar lz        = 4000.;    // m
735ccaff030SJeremy L Thompson   CeedScalar rc         = 1000.;    // m (Radius of bubble)
736ccaff030SJeremy L Thompson   PetscScalar resx      = 1000.;    // m (resolution in x)
737ccaff030SJeremy L Thompson   PetscScalar resy      = 1000.;    // m (resolution in y)
738ccaff030SJeremy L Thompson   PetscScalar resz      = 1000.;    // m (resolution in z)
739ccaff030SJeremy L Thompson   PetscInt outputfreq   = 10;       // -
740ccaff030SJeremy L Thompson   PetscInt contsteps    = 0;        // -
741ea6e0f84SLeila Ghaffari   PetscInt degreeVol    = 1;        // -
742ea6e0f84SLeila Ghaffari   PetscInt degreeSur    = 1;        // -
743ea6e0f84SLeila Ghaffari   PetscInt qextraVol    = 2;        // -
744ea6e0f84SLeila Ghaffari   PetscInt qextraSur    = 2;        // -
7450c6c0b13SLeila Ghaffari   DMBoundaryType periodicity[] = {DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,
7460c6c0b13SLeila Ghaffari                                   DM_BOUNDARY_NONE
7470c6c0b13SLeila Ghaffari                                  };
748ccaff030SJeremy L Thompson   PetscReal center[3], dc_axis[3] = {0, 0, 0};
749ccaff030SJeremy L Thompson 
750ccaff030SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
751ccaff030SJeremy L Thompson   if (ierr) return ierr;
752ccaff030SJeremy L Thompson 
753ccaff030SJeremy L Thompson   // Allocate PETSc context
754ccaff030SJeremy L Thompson   ierr = PetscCalloc1(1, &user); CHKERRQ(ierr);
755ccaff030SJeremy L Thompson   ierr = PetscMalloc1(1, &units); CHKERRQ(ierr);
756ccaff030SJeremy L Thompson 
757ccaff030SJeremy L Thompson   // Parse command line options
758ccaff030SJeremy L Thompson   comm = PETSC_COMM_WORLD;
759ccaff030SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
760ccaff030SJeremy L Thompson                            NULL); CHKERRQ(ierr);
761ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
762ccaff030SJeremy L Thompson                             NULL, ceedresource, ceedresource,
763ccaff030SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
7640c6c0b13SLeila Ghaffari   ierr = PetscOptionsBool("-test", "Run in test mode",
7650c6c0b13SLeila Ghaffari                           NULL, test=PETSC_FALSE, &test, NULL); CHKERRQ(ierr);
766ccaff030SJeremy L Thompson   problemChoice = NS_DENSITY_CURRENT;
767ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-problem", "Problem to solve", NULL,
768ccaff030SJeremy L Thompson                           problemTypes, (PetscEnum)problemChoice,
769ccaff030SJeremy L Thompson                           (PetscEnum *)&problemChoice, NULL); CHKERRQ(ierr);
770ccaff030SJeremy L Thompson   problem = &problemOptions[problemChoice];
771ccaff030SJeremy L Thompson   ierr = PetscOptionsEnum("-stab", "Stabilization method", NULL,
772ccaff030SJeremy L Thompson                           StabilizationTypes, (PetscEnum)(stab = STAB_NONE),
773ccaff030SJeremy L Thompson                           (PetscEnum *)&stab, NULL); CHKERRQ(ierr);
774ccaff030SJeremy L Thompson   ierr = PetscOptionsBool("-implicit", "Use implicit (IFunction) formulation",
775ccaff030SJeremy L Thompson                           NULL, implicit=PETSC_FALSE, &implicit, NULL);
776ccaff030SJeremy L Thompson   CHKERRQ(ierr);
777ccaff030SJeremy L Thompson   {
778ccaff030SJeremy L Thompson     PetscInt len;
779ccaff030SJeremy L Thompson     PetscBool flg;
780ccaff030SJeremy L Thompson     ierr = PetscOptionsIntArray("-bc_wall",
781ccaff030SJeremy L Thompson                                 "Use wall boundary conditions on this list of faces",
782ccaff030SJeremy L Thompson                                 NULL, bc.walls,
783ccaff030SJeremy L Thompson                                 (len = sizeof(bc.walls) / sizeof(bc.walls[0]),
784ccaff030SJeremy L Thompson                                  &len), &flg); CHKERRQ(ierr);
785ccaff030SJeremy L Thompson     if (flg) bc.nwall = len;
786ccaff030SJeremy L Thompson     for (PetscInt j=0; j<3; j++) {
787ccaff030SJeremy L Thompson       const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
788ccaff030SJeremy L Thompson       ierr = PetscOptionsIntArray(flags[j],
789ccaff030SJeremy L Thompson                                   "Use slip boundary conditions on this list of faces",
790ccaff030SJeremy L Thompson                                   NULL, bc.slips[j],
791ccaff030SJeremy L Thompson                                   (len = sizeof(bc.slips[j]) / sizeof(bc.slips[j][0]),
792ccaff030SJeremy L Thompson                                    &len), &flg);
793ccaff030SJeremy L Thompson       CHKERRQ(ierr);
7940c6c0b13SLeila Ghaffari       if (flg) bc.nslip[j] = len;
795ccaff030SJeremy L Thompson     }
796ccaff030SJeremy L Thompson   }
797cb3e2689Svaleriabarra   ierr = PetscOptionsInt("-viz_refine",
798cb3e2689Svaleriabarra                          "Regular refinement levels for visualization",
799cb3e2689Svaleriabarra                          NULL, viz_refine, &viz_refine, NULL);
800ccaff030SJeremy L Thompson   CHKERRQ(ierr);
801ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_meter", "1 meter in scaled length units",
802ccaff030SJeremy L Thompson                             NULL, meter, &meter, NULL); CHKERRQ(ierr);
803ccaff030SJeremy L Thompson   meter = fabs(meter);
804ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_second","1 second in scaled time units",
805ccaff030SJeremy L Thompson                             NULL, second, &second, NULL); CHKERRQ(ierr);
806ccaff030SJeremy L Thompson   second = fabs(second);
807ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_kilogram","1 kilogram in scaled mass units",
808ccaff030SJeremy L Thompson                             NULL, kilogram, &kilogram, NULL); CHKERRQ(ierr);
809ccaff030SJeremy L Thompson   kilogram = fabs(kilogram);
810ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-units_Kelvin",
811ccaff030SJeremy L Thompson                             "1 Kelvin in scaled temperature units",
812ccaff030SJeremy L Thompson                             NULL, Kelvin, &Kelvin, NULL); CHKERRQ(ierr);
813ccaff030SJeremy L Thompson   Kelvin = fabs(Kelvin);
814ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-theta0", "Reference potential temperature",
815ccaff030SJeremy L Thompson                             NULL, theta0, &theta0, NULL); CHKERRQ(ierr);
816ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-thetaC", "Perturbation of potential temperature",
817ccaff030SJeremy L Thompson                             NULL, thetaC, &thetaC, NULL); CHKERRQ(ierr);
818ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-P0", "Atmospheric pressure",
819ccaff030SJeremy L Thompson                             NULL, P0, &P0, NULL); CHKERRQ(ierr);
820ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-N", "Brunt-Vaisala frequency",
821ccaff030SJeremy L Thompson                             NULL, N, &N, NULL); CHKERRQ(ierr);
822ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cv", "Heat capacity at constant volume",
823ccaff030SJeremy L Thompson                             NULL, cv, &cv, NULL); CHKERRQ(ierr);
824ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-cp", "Heat capacity at constant pressure",
825ccaff030SJeremy L Thompson                             NULL, cp, &cp, NULL); CHKERRQ(ierr);
826ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-g", "Gravitational acceleration",
827ccaff030SJeremy L Thompson                             NULL, g, &g, NULL); CHKERRQ(ierr);
828ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lambda",
829ccaff030SJeremy L Thompson                             "Stokes hypothesis second viscosity coefficient",
830ccaff030SJeremy L Thompson                             NULL, lambda, &lambda, NULL); CHKERRQ(ierr);
831ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-mu", "Shear dynamic viscosity coefficient",
832ccaff030SJeremy L Thompson                             NULL, mu, &mu, NULL); CHKERRQ(ierr);
833ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-k", "Thermal conductivity",
834ccaff030SJeremy L Thompson                             NULL, k, &k, NULL); CHKERRQ(ierr);
835ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-CtauS",
836ccaff030SJeremy L Thompson                             "Scale coefficient for tau (nondimensional)",
837ccaff030SJeremy L Thompson                             NULL, CtauS, &CtauS, NULL); CHKERRQ(ierr);
838ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-strong_form",
839ccaff030SJeremy L Thompson                             "Strong (1) or weak/integrated by parts (0) advection residual",
840ccaff030SJeremy L Thompson                             NULL, strong_form, &strong_form, NULL);
841ccaff030SJeremy L Thompson   CHKERRQ(ierr);
842ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lx", "Length scale in x direction",
843ccaff030SJeremy L Thompson                             NULL, lx, &lx, NULL); CHKERRQ(ierr);
844ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-ly", "Length scale in y direction",
845ccaff030SJeremy L Thompson                             NULL, ly, &ly, NULL); CHKERRQ(ierr);
846ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-lz", "Length scale in z direction",
847ccaff030SJeremy L Thompson                             NULL, lz, &lz, NULL); CHKERRQ(ierr);
848ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-rc", "Characteristic radius of thermal bubble",
849ccaff030SJeremy L Thompson                             NULL, rc, &rc, NULL); CHKERRQ(ierr);
850ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resx","Target resolution in x",
851ccaff030SJeremy L Thompson                             NULL, resx, &resx, NULL); CHKERRQ(ierr);
852ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resy","Target resolution in y",
853ccaff030SJeremy L Thompson                             NULL, resy, &resy, NULL); CHKERRQ(ierr);
854ccaff030SJeremy L Thompson   ierr = PetscOptionsScalar("-resz","Target resolution in z",
855ccaff030SJeremy L Thompson                             NULL, resz, &resz, NULL); CHKERRQ(ierr);
856ccaff030SJeremy L Thompson   PetscInt n = problem->dim;
8570c6c0b13SLeila Ghaffari   ierr = PetscOptionsEnumArray("-periodicity", "Periodicity per direction",
8580c6c0b13SLeila Ghaffari                                NULL, DMBoundaryTypes, (PetscEnum *)periodicity,
8590c6c0b13SLeila Ghaffari                                &n, NULL); CHKERRQ(ierr);
8600c6c0b13SLeila Ghaffari   n = problem->dim;
861ccaff030SJeremy L Thompson   center[0] = 0.5 * lx;
862ccaff030SJeremy L Thompson   center[1] = 0.5 * ly;
863ccaff030SJeremy L Thompson   center[2] = 0.5 * lz;
864ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-center", "Location of bubble center",
865ccaff030SJeremy L Thompson                                NULL, center, &n, NULL); CHKERRQ(ierr);
866ccaff030SJeremy L Thompson   n = problem->dim;
867ccaff030SJeremy L Thompson   ierr = PetscOptionsRealArray("-dc_axis",
868ccaff030SJeremy L Thompson                                "Axis of density current cylindrical anomaly, or {0,0,0} for spherically symmetric",
869ccaff030SJeremy L Thompson                                NULL, dc_axis, &n, NULL); CHKERRQ(ierr);
870ccaff030SJeremy L Thompson   {
871ccaff030SJeremy L Thompson     PetscReal norm = PetscSqrtReal(PetscSqr(dc_axis[0]) +
872ccaff030SJeremy L Thompson                                    PetscSqr(dc_axis[1]) + PetscSqr(dc_axis[2]));
873ccaff030SJeremy L Thompson     if (norm > 0) {
874ccaff030SJeremy L Thompson       for (int i=0; i<3; i++) dc_axis[i] /= norm;
875ccaff030SJeremy L Thompson     }
876ccaff030SJeremy L Thompson   }
877ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-output_freq",
878ccaff030SJeremy L Thompson                          "Frequency of output, in number of steps",
879ccaff030SJeremy L Thompson                          NULL, outputfreq, &outputfreq, NULL); CHKERRQ(ierr);
880ccaff030SJeremy L Thompson   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
881ccaff030SJeremy L Thompson                          NULL, contsteps, &contsteps, NULL); CHKERRQ(ierr);
882ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-degree_volume", "Polynomial degree of finite elements for the Volume",
883ea6e0f84SLeila Ghaffari                          NULL, degreeVol, &degreeVol, NULL); CHKERRQ(ierr);
884ea6e0f84SLeila Ghaffari   ierr = PetscOptionsInt("-qextra_volume", "Number of extra quadrature points for the Volume",
885ea6e0f84SLeila Ghaffari                          NULL, qextraVol, &qextraVol, NULL); CHKERRQ(ierr);
8860c6c0b13SLeila Ghaffari   PetscStrncpy(user->outputfolder, ".", 2);
887ccaff030SJeremy L Thompson   ierr = PetscOptionsString("-of", "Output folder",
888ccaff030SJeremy L Thompson                             NULL, user->outputfolder, user->outputfolder,
889ccaff030SJeremy L Thompson                             sizeof(user->outputfolder), NULL); CHKERRQ(ierr);
890ccaff030SJeremy L Thompson   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
891ccaff030SJeremy L Thompson 
892ccaff030SJeremy L Thompson   // Define derived units
893ccaff030SJeremy L Thompson   Pascal = kilogram / (meter * PetscSqr(second));
894ccaff030SJeremy L Thompson   JperkgK =  PetscSqr(meter) / (PetscSqr(second) * Kelvin);
895ccaff030SJeremy L Thompson   mpersquareds = meter / PetscSqr(second);
896ccaff030SJeremy L Thompson   WpermK = kilogram * meter / (pow(second,3) * Kelvin);
897ccaff030SJeremy L Thompson   kgpercubicm = kilogram / pow(meter,3);
898ccaff030SJeremy L Thompson   kgpersquaredms = kilogram / (PetscSqr(meter) * second);
899ccaff030SJeremy L Thompson   Joulepercubicm = kilogram / (meter * PetscSqr(second));
900ccaff030SJeremy L Thompson 
901ccaff030SJeremy L Thompson   // Scale variables to desired units
902ccaff030SJeremy L Thompson   theta0 *= Kelvin;
903ccaff030SJeremy L Thompson   thetaC *= Kelvin;
904ccaff030SJeremy L Thompson   P0 *= Pascal;
905ccaff030SJeremy L Thompson   N *= (1./second);
906ccaff030SJeremy L Thompson   cv *= JperkgK;
907ccaff030SJeremy L Thompson   cp *= JperkgK;
908ccaff030SJeremy L Thompson   Rd = cp - cv;
909ccaff030SJeremy L Thompson   g *= mpersquareds;
910ccaff030SJeremy L Thompson   mu *= Pascal * second;
911ccaff030SJeremy L Thompson   k *= WpermK;
912ccaff030SJeremy L Thompson   lx = fabs(lx) * meter;
913ccaff030SJeremy L Thompson   ly = fabs(ly) * meter;
914ccaff030SJeremy L Thompson   lz = fabs(lz) * meter;
915ccaff030SJeremy L Thompson   rc = fabs(rc) * meter;
916ccaff030SJeremy L Thompson   resx = fabs(resx) * meter;
917ccaff030SJeremy L Thompson   resy = fabs(resy) * meter;
918ccaff030SJeremy L Thompson   resz = fabs(resz) * meter;
919ccaff030SJeremy L Thompson   for (int i=0; i<3; i++) center[i] *= meter;
920ccaff030SJeremy L Thompson 
921ccaff030SJeremy L Thompson   const CeedInt dim = problem->dim, ncompx = problem->dim,
922c96c872fSLeila Ghaffari                 qdatasizeVol = problem->qdatasizeVol,
923c96c872fSLeila Ghaffari                 qdatasizeSur = problem->qdatasizeSur;
924ccaff030SJeremy L Thompson   // Set up the libCEED context
925ccaff030SJeremy L Thompson   struct SetupContext_ ctxSetup =  {
926ccaff030SJeremy L Thompson     .theta0 = theta0,
927ccaff030SJeremy L Thompson     .thetaC = thetaC,
928ccaff030SJeremy L Thompson     .P0 = P0,
929ccaff030SJeremy L Thompson     .N = N,
930ccaff030SJeremy L Thompson     .cv = cv,
931ccaff030SJeremy L Thompson     .cp = cp,
932ccaff030SJeremy L Thompson     .Rd = Rd,
933ccaff030SJeremy L Thompson     .g = g,
934ccaff030SJeremy L Thompson     .rc = rc,
935ccaff030SJeremy L Thompson     .lx = lx,
936ccaff030SJeremy L Thompson     .ly = ly,
937ccaff030SJeremy L Thompson     .lz = lz,
9380c6c0b13SLeila Ghaffari     .periodicity0 = periodicity[0],
9390c6c0b13SLeila Ghaffari     .periodicity1 = periodicity[1],
9400c6c0b13SLeila Ghaffari     .periodicity2 = periodicity[2],
941ccaff030SJeremy L Thompson     .center[0] = center[0],
942ccaff030SJeremy L Thompson     .center[1] = center[1],
943ccaff030SJeremy L Thompson     .center[2] = center[2],
944ccaff030SJeremy L Thompson     .dc_axis[0] = dc_axis[0],
945ccaff030SJeremy L Thompson     .dc_axis[1] = dc_axis[1],
946ccaff030SJeremy L Thompson     .dc_axis[2] = dc_axis[2],
947ccaff030SJeremy L Thompson     .time = 0,
948ccaff030SJeremy L Thompson   };
949ccaff030SJeremy L Thompson 
950ccaff030SJeremy L Thompson   {
951ccaff030SJeremy L Thompson     const PetscReal scale[3] = {lx, ly, lz};
952ccaff030SJeremy L Thompson     ierr = DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, NULL, NULL, scale,
9530c6c0b13SLeila Ghaffari                                periodicity, PETSC_TRUE, &dm);
954ccaff030SJeremy L Thompson     CHKERRQ(ierr);
955ccaff030SJeremy L Thompson   }
9560c6c0b13SLeila Ghaffari   if (1) {
957ccaff030SJeremy L Thompson     DM               dmDist = NULL;
958ccaff030SJeremy L Thompson     PetscPartitioner part;
959ccaff030SJeremy L Thompson 
960ccaff030SJeremy L Thompson     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
961ccaff030SJeremy L Thompson     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
962ccaff030SJeremy L Thompson     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
963ccaff030SJeremy L Thompson     if (dmDist) {
964ccaff030SJeremy L Thompson       ierr = DMDestroy(&dm); CHKERRQ(ierr);
965ccaff030SJeremy L Thompson       dm  = dmDist;
966ccaff030SJeremy L Thompson     }
967ccaff030SJeremy L Thompson   }
968ccaff030SJeremy L Thompson   ierr = DMViewFromOptions(dm, NULL, "-dm_view"); CHKERRQ(ierr);
969ccaff030SJeremy L Thompson 
970ccaff030SJeremy L Thompson   ierr = DMLocalizeCoordinates(dm); CHKERRQ(ierr);
971ccaff030SJeremy L Thompson   ierr = DMSetFromOptions(dm); CHKERRQ(ierr);
972ea6e0f84SLeila Ghaffari   ierr = SetUpDM(dm, problem, degreeVol, &bc, &ctxSetup); CHKERRQ(ierr);
9730c6c0b13SLeila Ghaffari   if (!test) {
9740c6c0b13SLeila Ghaffari     ierr = PetscPrintf(PETSC_COMM_WORLD,
975ea6e0f84SLeila Ghaffari                        "Degree of Volumetric FEM Space: %D\n",
976ea6e0f84SLeila Ghaffari                        degreeVol); CHKERRQ(ierr);
9770c6c0b13SLeila Ghaffari   }
978ccaff030SJeremy L Thompson   dmviz = NULL;
979ccaff030SJeremy L Thompson   interpviz = NULL;
980ccaff030SJeremy L Thompson   if (viz_refine) {
981ff6701fcSJed Brown     DM dmhierarchy[viz_refine+1];
982ff6701fcSJed Brown 
983ccaff030SJeremy L Thompson     ierr = DMPlexSetRefinementUniform(dm, PETSC_TRUE); CHKERRQ(ierr);
984ff6701fcSJed Brown     dmhierarchy[0] = dm;
985ea6e0f84SLeila Ghaffari     for (PetscInt i = 0, d = degreeVol; i < viz_refine; i++) {
986ff6701fcSJed Brown       Mat interp_next;
987ff6701fcSJed Brown 
988ff6701fcSJed Brown       ierr = DMRefine(dmhierarchy[i], MPI_COMM_NULL, &dmhierarchy[i+1]);
989ccaff030SJeremy L Thompson       CHKERRQ(ierr);
990ff6701fcSJed Brown       ierr = DMSetCoarseDM(dmhierarchy[i+1], dmhierarchy[i]); CHKERRQ(ierr);
991ff6701fcSJed Brown       d = (d + 1) / 2;
992ff6701fcSJed Brown       if (i + 1 == viz_refine) d = 1;
993ff6701fcSJed Brown       ierr = SetUpDM(dmhierarchy[i+1], problem, d, &bc, &ctxSetup); CHKERRQ(ierr);
994ff6701fcSJed Brown       ierr = DMCreateInterpolation(dmhierarchy[i], dmhierarchy[i+1],
995ff6701fcSJed Brown                                    &interp_next, NULL); CHKERRQ(ierr);
996ff6701fcSJed Brown       if (!i) interpviz = interp_next;
997ff6701fcSJed Brown       else {
998ff6701fcSJed Brown         Mat C;
999ff6701fcSJed Brown         ierr = MatMatMult(interp_next, interpviz, MAT_INITIAL_MATRIX,
1000ff6701fcSJed Brown                           PETSC_DECIDE, &C); CHKERRQ(ierr);
1001ff6701fcSJed Brown         ierr = MatDestroy(&interp_next); CHKERRQ(ierr);
1002ff6701fcSJed Brown         ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1003ff6701fcSJed Brown         interpviz = C;
1004ff6701fcSJed Brown       }
1005ff6701fcSJed Brown     }
1006cb3e2689Svaleriabarra     for (PetscInt i=1; i<viz_refine; i++) {
1007ff6701fcSJed Brown       ierr = DMDestroy(&dmhierarchy[i]); CHKERRQ(ierr);
1008cb3e2689Svaleriabarra     }
1009ff6701fcSJed Brown     dmviz = dmhierarchy[viz_refine];
1010ccaff030SJeremy L Thompson   }
1011ccaff030SJeremy L Thompson   ierr = DMCreateGlobalVector(dm, &Q); CHKERRQ(ierr);
1012ccaff030SJeremy L Thompson   ierr = DMGetLocalVector(dm, &Qloc); CHKERRQ(ierr);
1013ccaff030SJeremy L Thompson   ierr = VecGetSize(Qloc, &lnodes); CHKERRQ(ierr);
1014ccaff030SJeremy L Thompson   lnodes /= ncompq;
1015ccaff030SJeremy L Thompson 
10160c6c0b13SLeila Ghaffari   {
10170c6c0b13SLeila Ghaffari     // Print grid information
1018ccaff030SJeremy L Thompson     CeedInt gdofs, odofs;
1019ccaff030SJeremy L Thompson     int comm_size;
1020ccaff030SJeremy L Thompson     char box_faces_str[PETSC_MAX_PATH_LEN] = "NONE";
1021ccaff030SJeremy L Thompson     ierr = VecGetSize(Q, &gdofs); CHKERRQ(ierr);
1022ccaff030SJeremy L Thompson     ierr = VecGetLocalSize(Q, &odofs); CHKERRQ(ierr);
1023ccaff030SJeremy L Thompson     ierr = MPI_Comm_size(comm, &comm_size); CHKERRQ(ierr);
1024ccaff030SJeremy L Thompson     ierr = PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str,
1025ccaff030SJeremy L Thompson                                  sizeof(box_faces_str), NULL); CHKERRQ(ierr);
10260c6c0b13SLeila Ghaffari     if (!test) {
10270c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Global FEM dofs: %D (%D owned) on %d ranks\n",
10280c6c0b13SLeila Ghaffari                          gdofs, odofs, comm_size); CHKERRQ(ierr);
10290c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "Local FEM nodes: %D\n", lnodes); CHKERRQ(ierr);
10300c6c0b13SLeila Ghaffari       ierr = PetscPrintf(comm, "dm_plex_box_faces: %s\n", box_faces_str);
1031ccaff030SJeremy L Thompson       CHKERRQ(ierr);
1032ccaff030SJeremy L Thompson     }
1033ccaff030SJeremy L Thompson 
10340c6c0b13SLeila Ghaffari   }
10350c6c0b13SLeila Ghaffari 
1036ccaff030SJeremy L Thompson   // Set up global mass vector
1037ccaff030SJeremy L Thompson   ierr = VecDuplicate(Q,&user->M); CHKERRQ(ierr);
1038ccaff030SJeremy L Thompson 
10390c6c0b13SLeila Ghaffari   // Set up CEED
1040ccaff030SJeremy L Thompson   // CEED Bases
1041ccaff030SJeremy L Thompson   CeedInit(ceedresource, &ceed);
1042c96c872fSLeila Ghaffari   // Bases for the Volume
1043ea6e0f84SLeila Ghaffari   numP_Vol = degreeVol + 1;
1044ea6e0f84SLeila Ghaffari   numQ_Vol = numP_Vol + qextraVol;
1045ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompq, numP_Vol, numQ_Vol, CEED_GAUSS,
1046ea6e0f84SLeila Ghaffari                                   &basisqVol);
1047ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numQ_Vol, CEED_GAUSS,
1048ea6e0f84SLeila Ghaffari                                   &basisxVol);
1049ea6e0f84SLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim, ncompx, 2, numP_Vol,
1050ea6e0f84SLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcVol);
1051c96c872fSLeila Ghaffari   // Bases for the Surface
1052c96c872fSLeila Ghaffari   CeedInt height = 1;
1053c96c872fSLeila Ghaffari   numP_Sur = degreeSur + 1;
1054c96c872fSLeila Ghaffari   numQ_Sur = numP_Sur + qextraSur;
1055c96c872fSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim - height, ncompq, numP_Sur, numQ_Sur, CEED_GAUSS,
1056c96c872fSLeila Ghaffari                                   &basisqSur);
1057c96c872fSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim - height, ncompx, 2, numQ_Sur, CEED_GAUSS,
1058c96c872fSLeila Ghaffari                                   &basisxSur);
1059c96c872fSLeila Ghaffari   CeedBasisCreateTensorH1Lagrange(ceed, dim - height, ncompx, 2, numP_Sur,
1060c96c872fSLeila Ghaffari                                   CEED_GAUSS_LOBATTO, &basisxcSur);
1061ccaff030SJeremy L Thompson 
1062ccaff030SJeremy L Thompson   ierr = DMGetCoordinateDM(dm, &dmcoord); CHKERRQ(ierr);
1063ccaff030SJeremy L Thompson   ierr = DMPlexSetClosurePermutationTensor(dmcoord, PETSC_DETERMINE, NULL);
1064ccaff030SJeremy L Thompson   CHKERRQ(ierr);
1065ccaff030SJeremy L Thompson 
1066ccaff030SJeremy L Thompson   // CEED Restrictions
1067c96c872fSLeila Ghaffari   // Restrictions on the Volume
1068c96c872fSLeila Ghaffari   ierr = GetRestriction(ceed, dm, ncompx, dim, 0, 0, 0, numP_Vol, numQ_Vol, qdatasizeVol,
1069c96c872fSLeila Ghaffari     &restrictqVol, &restrictxVol, &restrictqdiVol); CHKERRQ(ierr);
1070c96c872fSLeila Ghaffari   // (Rough draft) Restriction for one face ----> Should be done for all Neumann faces
1071c96c872fSLeila Ghaffari   ierr = GetRestriction(ceed, dm, ncompx, dim - height, height, "Face Sets", 6, numP_Sur,
1072c96c872fSLeila Ghaffari     numQ_Sur, qdatasizeSur, &restrictqSur, &restrictxSur, &restrictqdiSur); CHKERRQ(ierr);
1073ccaff030SJeremy L Thompson 
1074ccaff030SJeremy L Thompson   ierr = DMGetCoordinatesLocal(dm, &Xloc); CHKERRQ(ierr);
1075ccaff030SJeremy L Thompson   ierr = CreateVectorFromPetscVec(ceed, Xloc, &xcorners); CHKERRQ(ierr);
1076ccaff030SJeremy L Thompson 
1077ccaff030SJeremy L Thompson   // Create the CEED vectors that will be needed in setup
1078c96c872fSLeila Ghaffari   CeedInt NqptsVol, NqptsSur;
1079c96c872fSLeila Ghaffari   // Volume
1080c96c872fSLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqVol, &NqptsVol);
1081c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictqVol, &localNelemVol);
1082c96c872fSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeVol*localNelemVol*NqptsVol, &qdataVol);
1083c96c872fSLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &q0ceedVol, NULL);
1084c96c872fSLeila Ghaffari   // Surface
1085c96c872fSLeila Ghaffari   CeedBasisGetNumQuadraturePoints(basisqSur, &NqptsSur);
1086c96c872fSLeila Ghaffari   CeedElemRestrictionGetNumElements(restrictqSur, &localNelemSur);
1087c96c872fSLeila Ghaffari   CeedVectorCreate(ceed, qdatasizeSur*localNelemSur*NqptsSur, &qdataSur);
1088c96c872fSLeila Ghaffari   //CeedElemRestrictionCreateVector(restrictqSur, &q0ceedSur, NULL);
1089ccaff030SJeremy L Thompson 
1090ccaff030SJeremy L Thompson   // Create the Q-function that builds the quadrature data for the NS operator
1091ea6e0f84SLeila Ghaffari   CeedQFunctionCreateInterior(ceed, 1, problem->setupVol, problem->setupVol_loc,
1092ea6e0f84SLeila Ghaffari                               &qf_setupVol);
1093ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "dx", ncompx*dim, CEED_EVAL_GRAD);
1094ea6e0f84SLeila Ghaffari   CeedQFunctionAddInput(qf_setupVol, "weight", 1, CEED_EVAL_WEIGHT);
1095ea6e0f84SLeila Ghaffari   CeedQFunctionAddOutput(qf_setupVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1096ccaff030SJeremy L Thompson 
1097ccaff030SJeremy L Thompson   // Create the Q-function that sets the ICs of the operator
1098ccaff030SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics, problem->ics_loc, &qf_ics);
1099ccaff030SJeremy L Thompson   CeedQFunctionAddInput(qf_ics, "x", ncompx, CEED_EVAL_INTERP);
1100ccaff030SJeremy L Thompson   CeedQFunctionAddOutput(qf_ics, "q0", ncompq, CEED_EVAL_NONE);
1101ccaff030SJeremy L Thompson 
1102ea6e0f84SLeila Ghaffari   qf_rhsVol = NULL;
1103ea6e0f84SLeila Ghaffari   if (problem->applyVol_rhs) { // Create the Q-function that defines the action of the RHS operator
1104ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_rhs,
1105ea6e0f84SLeila Ghaffari                                 problem->applyVol_rhs_loc, &qf_rhsVol);
1106ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "q", ncompq, CEED_EVAL_INTERP);
1107ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1108ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1109ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_rhsVol, "x", ncompx, CEED_EVAL_INTERP);
1110ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "v", ncompq, CEED_EVAL_INTERP);
1111ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_rhsVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1112ccaff030SJeremy L Thompson   }
1113ccaff030SJeremy L Thompson 
1114ea6e0f84SLeila Ghaffari   qf_ifunctionVol = NULL;
1115ea6e0f84SLeila Ghaffari   if (problem->applyVol_ifunction) { // Create the Q-function that defines the action of the IFunction
1116ea6e0f84SLeila Ghaffari     CeedQFunctionCreateInterior(ceed, 1, problem->applyVol_ifunction,
1117ea6e0f84SLeila Ghaffari                                 problem->applyVol_ifunction_loc, &qf_ifunctionVol);
1118ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "q", ncompq, CEED_EVAL_INTERP);
1119ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "dq", ncompq*dim, CEED_EVAL_GRAD);
1120ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdot", ncompq, CEED_EVAL_INTERP);
1121ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "qdataVol", qdatasizeVol, CEED_EVAL_NONE);
1122ea6e0f84SLeila Ghaffari     CeedQFunctionAddInput(qf_ifunctionVol, "x", ncompx, CEED_EVAL_INTERP);
1123ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "v", ncompq, CEED_EVAL_INTERP);
1124ea6e0f84SLeila Ghaffari     CeedQFunctionAddOutput(qf_ifunctionVol, "dv", ncompq*dim, CEED_EVAL_GRAD);
1125ccaff030SJeremy L Thompson   }
1126ccaff030SJeremy L Thompson 
1127ccaff030SJeremy L Thompson   // Create the operator that builds the quadrature data for the NS operator
1128ea6e0f84SLeila Ghaffari   CeedOperatorCreate(ceed, qf_setupVol, NULL, NULL, &op_setupVol);
1129ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "dx", restrictxVol, basisxVol, CEED_VECTOR_ACTIVE);
1130ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "weight", CEED_ELEMRESTRICTION_NONE,
1131ea6e0f84SLeila Ghaffari                        basisxVol, CEED_VECTOR_NONE);
1132ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_setupVol, "qdataVol", restrictqdiVol,
1133ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1134ccaff030SJeremy L Thompson 
1135ccaff030SJeremy L Thompson   // Create the operator that sets the ICs
1136ccaff030SJeremy L Thompson   CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics);
1137ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "x", restrictxVol, basisxcVol, CEED_VECTOR_ACTIVE);
1138ea6e0f84SLeila Ghaffari   CeedOperatorSetField(op_ics, "q0", restrictqVol,
1139ccaff030SJeremy L Thompson                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1140ccaff030SJeremy L Thompson 
1141ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qceed, NULL);
1142ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->qdotceed, NULL);
1143ea6e0f84SLeila Ghaffari   CeedElemRestrictionCreateVector(restrictqVol, &user->gceed, NULL);
1144ccaff030SJeremy L Thompson 
1145ea6e0f84SLeila Ghaffari   if (qf_rhsVol) { // Create the RHS physics operator
1146ccaff030SJeremy L Thompson     CeedOperator op;
1147ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_rhsVol, NULL, NULL, &op);
1148ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1149ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1150ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdataVol", restrictqdiVol,
1151ea6e0f84SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataVol);
1152ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1153ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1154ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1155ccaff030SJeremy L Thompson     user->op_rhs = op;
1156ccaff030SJeremy L Thompson   }
1157ccaff030SJeremy L Thompson 
1158ea6e0f84SLeila Ghaffari   if (qf_ifunctionVol) { // Create the IFunction operator
1159ccaff030SJeremy L Thompson     CeedOperator op;
1160ea6e0f84SLeila Ghaffari     CeedOperatorCreate(ceed, qf_ifunctionVol, NULL, NULL, &op);
1161ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "q", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1162ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dq", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1163ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdot", restrictqVol, basisqVol, user->qdotceed);
1164ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "qdataVol", restrictqdiVol,
1165ea6e0f84SLeila Ghaffari                          CEED_BASIS_COLLOCATED, qdataVol);
1166ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "x", restrictxVol, basisxVol, xcorners);
1167ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "v", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1168ea6e0f84SLeila Ghaffari     CeedOperatorSetField(op, "dv", restrictqVol, basisqVol, CEED_VECTOR_ACTIVE);
1169ccaff030SJeremy L Thompson     user->op_ifunction = op;
1170ccaff030SJeremy L Thompson   }
1171ccaff030SJeremy L Thompson 
1172ccaff030SJeremy L Thompson   CeedQFunctionSetContext(qf_ics, &ctxSetup, sizeof ctxSetup);
1173ccaff030SJeremy L Thompson   CeedScalar ctxNS[8] = {lambda, mu, k, cv, cp, g, Rd};
1174ccaff030SJeremy L Thompson   struct Advection2dContext_ ctxAdvection2d = {
1175ccaff030SJeremy L Thompson     .CtauS = CtauS,
1176ccaff030SJeremy L Thompson     .strong_form = strong_form,
1177ccaff030SJeremy L Thompson     .stabilization = stab,
1178ccaff030SJeremy L Thompson   };
1179ccaff030SJeremy L Thompson   switch (problemChoice) {
1180ccaff030SJeremy L Thompson   case NS_DENSITY_CURRENT:
1181ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxNS, sizeof ctxNS);
1182ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxNS,
1183ccaff030SJeremy L Thompson           sizeof ctxNS);
1184ccaff030SJeremy L Thompson     break;
1185ccaff030SJeremy L Thompson   case NS_ADVECTION:
1186ccaff030SJeremy L Thompson   case NS_ADVECTION2D:
1187ea6e0f84SLeila Ghaffari     if (qf_rhsVol) CeedQFunctionSetContext(qf_rhsVol, &ctxAdvection2d,
1188ccaff030SJeremy L Thompson                                           sizeof ctxAdvection2d);
1189ea6e0f84SLeila Ghaffari     if (qf_ifunctionVol) CeedQFunctionSetContext(qf_ifunctionVol, &ctxAdvection2d,
1190ccaff030SJeremy L Thompson           sizeof ctxAdvection2d);
1191ccaff030SJeremy L Thompson   }
1192ccaff030SJeremy L Thompson 
1193ccaff030SJeremy L Thompson   // Set up PETSc context
1194ccaff030SJeremy L Thompson   // Set up units structure
1195ccaff030SJeremy L Thompson   units->meter = meter;
1196ccaff030SJeremy L Thompson   units->kilogram = kilogram;
1197ccaff030SJeremy L Thompson   units->second = second;
1198ccaff030SJeremy L Thompson   units->Kelvin = Kelvin;
1199ccaff030SJeremy L Thompson   units->Pascal = Pascal;
1200ccaff030SJeremy L Thompson   units->JperkgK = JperkgK;
1201ccaff030SJeremy L Thompson   units->mpersquareds = mpersquareds;
1202ccaff030SJeremy L Thompson   units->WpermK = WpermK;
1203ccaff030SJeremy L Thompson   units->kgpercubicm = kgpercubicm;
1204ccaff030SJeremy L Thompson   units->kgpersquaredms = kgpersquaredms;
1205ccaff030SJeremy L Thompson   units->Joulepercubicm = Joulepercubicm;
1206ccaff030SJeremy L Thompson 
1207ccaff030SJeremy L Thompson   // Set up user structure
1208ccaff030SJeremy L Thompson   user->comm = comm;
1209ccaff030SJeremy L Thompson   user->outputfreq = outputfreq;
1210ccaff030SJeremy L Thompson   user->contsteps = contsteps;
1211ccaff030SJeremy L Thompson   user->units = units;
1212ccaff030SJeremy L Thompson   user->dm = dm;
1213ccaff030SJeremy L Thompson   user->dmviz = dmviz;
1214ccaff030SJeremy L Thompson   user->interpviz = interpviz;
1215ccaff030SJeremy L Thompson   user->ceed = ceed;
1216ccaff030SJeremy L Thompson 
1217ea6e0f84SLeila Ghaffari   // Calculate qdataVol and ICs
1218ccaff030SJeremy L Thompson   // Set up state global and local vectors
1219ccaff030SJeremy L Thompson   ierr = VecZeroEntries(Q); CHKERRQ(ierr);
1220ccaff030SJeremy L Thompson 
1221c96c872fSLeila Ghaffari   ierr = VectorPlacePetscVec(q0ceedVol, Qloc); CHKERRQ(ierr);
1222ccaff030SJeremy L Thompson 
1223ccaff030SJeremy L Thompson   // Apply Setup Ceed Operators
1224ccaff030SJeremy L Thompson   ierr = VectorPlacePetscVec(xcorners, Xloc); CHKERRQ(ierr);
1225ea6e0f84SLeila Ghaffari   CeedOperatorApply(op_setupVol, xcorners, qdataVol, CEED_REQUEST_IMMEDIATE);
1226ea6e0f84SLeila Ghaffari   ierr = ComputeLumpedMassMatrix(ceed, dm, restrictqVol, basisqVol, restrictqdiVol, qdataVol,
1227ccaff030SJeremy L Thompson                                  user->M); CHKERRQ(ierr);
1228ccaff030SJeremy L Thompson 
1229c96c872fSLeila Ghaffari   ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceedVol, dm, Qloc, Q, restrictqVol,
12300c6c0b13SLeila Ghaffari                                &ctxSetup, 0.0);
1231ccaff030SJeremy L Thompson   if (1) { // Record boundary values from initial condition and override DMPlexInsertBoundaryValues()
1232ccaff030SJeremy L Thompson     // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow.  If we
1233ccaff030SJeremy L Thompson     // disable this, we should still get the same results due to the problem->bc function, but with potentially much
1234ccaff030SJeremy L Thompson     // slower execution.
1235ccaff030SJeremy L Thompson     Vec Qbc;
1236ccaff030SJeremy L Thompson     ierr = DMGetNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1237ccaff030SJeremy L Thompson     ierr = VecCopy(Qloc, Qbc); CHKERRQ(ierr);
1238ccaff030SJeremy L Thompson     ierr = VecZeroEntries(Qloc); CHKERRQ(ierr);
1239ccaff030SJeremy L Thompson     ierr = DMGlobalToLocal(dm, Q, INSERT_VALUES, Qloc); CHKERRQ(ierr);
1240ccaff030SJeremy L Thompson     ierr = VecAXPY(Qbc, -1., Qloc); CHKERRQ(ierr);
1241ccaff030SJeremy L Thompson     ierr = DMRestoreNamedLocalVector(dm, "Qbc", &Qbc); CHKERRQ(ierr);
1242ccaff030SJeremy L Thompson     ierr = PetscObjectComposeFunction((PetscObject)dm,
12430c6c0b13SLeila Ghaffari                                       "DMPlexInsertBoundaryValues_C",DMPlexInsertBoundaryValues_NS); CHKERRQ(ierr);
1244ccaff030SJeremy L Thompson   }
1245ccaff030SJeremy L Thompson 
1246ccaff030SJeremy L Thompson   MPI_Comm_rank(comm, &rank);
1247ccaff030SJeremy L Thompson   if (!rank) {ierr = PetscMkdir(user->outputfolder); CHKERRQ(ierr);}
1248ccaff030SJeremy L Thompson   // Gather initial Q values
1249ccaff030SJeremy L Thompson   // In case of continuation of simulation, set up initial values from binary file
1250ccaff030SJeremy L Thompson   if (contsteps) { // continue from existent solution
1251ccaff030SJeremy L Thompson     PetscViewer viewer;
1252ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1253ccaff030SJeremy L Thompson     // Read input
1254ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-solution.bin",
1255ccaff030SJeremy L Thompson                          user->outputfolder);
1256ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1257ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1258ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1259ccaff030SJeremy L Thompson     ierr = VecLoad(Q, viewer); CHKERRQ(ierr);
1260ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
12610c6c0b13SLeila Ghaffari   } else {
12620c6c0b13SLeila Ghaffari     //ierr = DMLocalToGlobal(dm, Qloc, INSERT_VALUES, Q);CHKERRQ(ierr);
1263ccaff030SJeremy L Thompson   }
1264ccaff030SJeremy L Thompson   ierr = DMRestoreLocalVector(dm, &Qloc); CHKERRQ(ierr);
1265ccaff030SJeremy L Thompson 
1266ccaff030SJeremy L Thompson   // Create and setup TS
1267ccaff030SJeremy L Thompson   ierr = TSCreate(comm, &ts); CHKERRQ(ierr);
1268ccaff030SJeremy L Thompson   ierr = TSSetDM(ts, dm); CHKERRQ(ierr);
1269ccaff030SJeremy L Thompson   if (implicit) {
1270ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSBDF); CHKERRQ(ierr);
1271ccaff030SJeremy L Thompson     if (user->op_ifunction) {
1272ccaff030SJeremy L Thompson       ierr = TSSetIFunction(ts, NULL, IFunction_NS, &user); CHKERRQ(ierr);
1273ccaff030SJeremy L Thompson     } else {                    // Implicit integrators can fall back to using an RHSFunction
1274ccaff030SJeremy L Thompson       ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1275ccaff030SJeremy L Thompson     }
1276ccaff030SJeremy L Thompson   } else {
1277ccaff030SJeremy L Thompson     if (!user->op_rhs) SETERRQ(comm,PETSC_ERR_ARG_NULL,
1278ccaff030SJeremy L Thompson                                  "Problem does not provide RHSFunction");
1279ccaff030SJeremy L Thompson     ierr = TSSetType(ts, TSRK); CHKERRQ(ierr);
1280ccaff030SJeremy L Thompson     ierr = TSRKSetType(ts, TSRK5F); CHKERRQ(ierr);
1281ccaff030SJeremy L Thompson     ierr = TSSetRHSFunction(ts, NULL, RHS_NS, &user); CHKERRQ(ierr);
1282ccaff030SJeremy L Thompson   }
1283ccaff030SJeremy L Thompson   ierr = TSSetMaxTime(ts, 500. * units->second); CHKERRQ(ierr);
1284ccaff030SJeremy L Thompson   ierr = TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER); CHKERRQ(ierr);
1285ccaff030SJeremy L Thompson   ierr = TSSetTimeStep(ts, 1.e-2 * units->second); CHKERRQ(ierr);
12860c6c0b13SLeila Ghaffari   if (test) {ierr = TSSetMaxSteps(ts, 1); CHKERRQ(ierr);}
1287ccaff030SJeremy L Thompson   ierr = TSGetAdapt(ts, &adapt); CHKERRQ(ierr);
1288ccaff030SJeremy L Thompson   ierr = TSAdaptSetStepLimits(adapt,
1289ccaff030SJeremy L Thompson                               1.e-12 * units->second,
1290ccaff030SJeremy L Thompson                               1.e2 * units->second); CHKERRQ(ierr);
1291ccaff030SJeremy L Thompson   ierr = TSSetFromOptions(ts); CHKERRQ(ierr);
1292ccaff030SJeremy L Thompson   if (!contsteps) { // print initial condition
12930c6c0b13SLeila Ghaffari     if (!test) {
1294ccaff030SJeremy L Thompson       ierr = TSMonitor_NS(ts, 0, 0., Q, user); CHKERRQ(ierr);
1295ccaff030SJeremy L Thompson     }
1296ccaff030SJeremy L Thompson   } else { // continue from time of last output
1297ccaff030SJeremy L Thompson     PetscReal time;
1298ccaff030SJeremy L Thompson     PetscInt count;
1299ccaff030SJeremy L Thompson     PetscViewer viewer;
1300ccaff030SJeremy L Thompson     char filepath[PETSC_MAX_PATH_LEN];
1301ccaff030SJeremy L Thompson     ierr = PetscSNPrintf(filepath, sizeof filepath, "%s/ns-time.bin",
1302ccaff030SJeremy L Thompson                          user->outputfolder); CHKERRQ(ierr);
1303ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryOpen(comm, filepath, FILE_MODE_READ, &viewer);
1304ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1305ccaff030SJeremy L Thompson     ierr = PetscViewerBinaryRead(viewer, &time, 1, &count, PETSC_REAL);
1306ccaff030SJeremy L Thompson     CHKERRQ(ierr);
1307ccaff030SJeremy L Thompson     ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr);
1308ccaff030SJeremy L Thompson     ierr = TSSetTime(ts, time * user->units->second); CHKERRQ(ierr);
1309ccaff030SJeremy L Thompson   }
13100c6c0b13SLeila Ghaffari   if (!test) {
1311ccaff030SJeremy L Thompson     ierr = TSMonitorSet(ts, TSMonitor_NS, user, NULL); CHKERRQ(ierr);
1312ccaff030SJeremy L Thompson   }
1313ccaff030SJeremy L Thompson 
1314ccaff030SJeremy L Thompson   // Solve
1315ccaff030SJeremy L Thompson   start = MPI_Wtime();
1316ccaff030SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ts); CHKERRQ(ierr);
1317ccaff030SJeremy L Thompson   ierr = TSSolve(ts, Q); CHKERRQ(ierr);
1318ccaff030SJeremy L Thompson   cpu_time_used = MPI_Wtime() - start;
1319ccaff030SJeremy L Thompson   ierr = TSGetSolveTime(ts,&ftime); CHKERRQ(ierr);
1320ccaff030SJeremy L Thompson   ierr = MPI_Allreduce(MPI_IN_PLACE, &cpu_time_used, 1, MPI_DOUBLE, MPI_MIN,
1321ccaff030SJeremy L Thompson                        comm); CHKERRQ(ierr);
13220c6c0b13SLeila Ghaffari   if (!test) {
1323ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
13240c6c0b13SLeila Ghaffari                        "Time taken for solution: %g\n",
1325ccaff030SJeremy L Thompson                        (double)cpu_time_used); CHKERRQ(ierr);
1326ccaff030SJeremy L Thompson   }
1327ccaff030SJeremy L Thompson 
1328ccaff030SJeremy L Thompson   // Get error
13290c6c0b13SLeila Ghaffari   if (problem->non_zero_time && !test) {
1330ccaff030SJeremy L Thompson     Vec Qexact, Qexactloc;
1331ccaff030SJeremy L Thompson     PetscReal norm;
1332ccaff030SJeremy L Thompson     ierr = DMCreateGlobalVector(dm, &Qexact); CHKERRQ(ierr);
1333ccaff030SJeremy L Thompson     ierr = DMGetLocalVector(dm, &Qexactloc); CHKERRQ(ierr);
1334ccaff030SJeremy L Thompson     ierr = VecGetSize(Qexactloc, &lnodes); CHKERRQ(ierr);
1335ccaff030SJeremy L Thompson 
1336c96c872fSLeila Ghaffari     ierr = ICs_PetscMultiplicity(op_ics, xcorners, q0ceedVol, dm, Qexactloc, Qexact,
1337ea6e0f84SLeila Ghaffari                                  restrictqVol, &ctxSetup, ftime); CHKERRQ(ierr);
1338ccaff030SJeremy L Thompson 
1339ccaff030SJeremy L Thompson     ierr = VecAXPY(Q, -1.0, Qexact);  CHKERRQ(ierr);
1340ccaff030SJeremy L Thompson     ierr = VecNorm(Q, NORM_MAX, &norm); CHKERRQ(ierr);
1341c96c872fSLeila Ghaffari     CeedVectorDestroy(&q0ceedVol);
1342ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1343ccaff030SJeremy L Thompson                        "Max Error: %g\n",
1344ccaff030SJeremy L Thompson                        (double)norm); CHKERRQ(ierr);
1345ccaff030SJeremy L Thompson   }
1346ccaff030SJeremy L Thompson 
1347ccaff030SJeremy L Thompson   // Output Statistics
1348ccaff030SJeremy L Thompson   ierr = TSGetStepNumber(ts,&steps); CHKERRQ(ierr);
13490c6c0b13SLeila Ghaffari   if (!test) {
1350ccaff030SJeremy L Thompson     ierr = PetscPrintf(PETSC_COMM_WORLD,
1351ccaff030SJeremy L Thompson                        "Time integrator took %D time steps to reach final time %g\n",
1352ccaff030SJeremy L Thompson                        steps,(double)ftime); CHKERRQ(ierr);
1353ccaff030SJeremy L Thompson   }
13549cf88b28Svaleriabarra 
1355ccaff030SJeremy L Thompson   // Clean up libCEED
1356ea6e0f84SLeila Ghaffari   CeedVectorDestroy(&qdataVol);
1357ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qceed);
1358ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->qdotceed);
1359ccaff030SJeremy L Thompson   CeedVectorDestroy(&user->gceed);
1360ccaff030SJeremy L Thompson   CeedVectorDestroy(&xcorners);
1361ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisqVol);
1362ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxVol);
1363ea6e0f84SLeila Ghaffari   CeedBasisDestroy(&basisxcVol);
1364ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqVol);
1365ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictxVol);
1366ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictqdiVol);
1367ea6e0f84SLeila Ghaffari   CeedElemRestrictionDestroy(&restrictxcoordVol);
1368ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_setupVol);
1369ccaff030SJeremy L Thompson   CeedQFunctionDestroy(&qf_ics);
1370ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_rhsVol);
1371ea6e0f84SLeila Ghaffari   CeedQFunctionDestroy(&qf_ifunctionVol);
1372ea6e0f84SLeila Ghaffari   CeedOperatorDestroy(&op_setupVol);
1373ccaff030SJeremy L Thompson   CeedOperatorDestroy(&op_ics);
1374ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_rhs);
1375ccaff030SJeremy L Thompson   CeedOperatorDestroy(&user->op_ifunction);
1376ccaff030SJeremy L Thompson   CeedDestroy(&ceed);
1377ccaff030SJeremy L Thompson 
1378ccaff030SJeremy L Thompson   // Clean up PETSc
1379ccaff030SJeremy L Thompson   ierr = VecDestroy(&Q); CHKERRQ(ierr);
1380ccaff030SJeremy L Thompson   ierr = VecDestroy(&user->M); CHKERRQ(ierr);
1381ccaff030SJeremy L Thompson   ierr = MatDestroy(&interpviz); CHKERRQ(ierr);
1382ccaff030SJeremy L Thompson   ierr = DMDestroy(&dmviz); CHKERRQ(ierr);
1383ccaff030SJeremy L Thompson   ierr = TSDestroy(&ts); CHKERRQ(ierr);
1384ccaff030SJeremy L Thompson   ierr = DMDestroy(&dm); CHKERRQ(ierr);
1385ccaff030SJeremy L Thompson   ierr = PetscFree(units); CHKERRQ(ierr);
1386ccaff030SJeremy L Thompson   ierr = PetscFree(user); CHKERRQ(ierr);
1387ccaff030SJeremy L Thompson   return PetscFinalize();
1388ccaff030SJeremy L Thompson }
1389ccaff030SJeremy L Thompson 
1390