xref: /libCEED/backends/opt/ceed-opt-operator.c (revision bd4df46207a63b48a485389a63848eb7a2477f6b)
1 // Copyright (c) 2017-2018, Lawrence Livermore National Security, LLC.
2 // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707.
3 // All Rights reserved. See files LICENSE and NOTICE for details.
4 //
5 // This file is part of CEED, a collection of benchmarks, miniapps, software
6 // libraries and APIs for efficient high-order finite element and spectral
7 // element discretizations for exascale applications. For more information and
8 // source code availability see http://github.com/ceed.
9 //
10 // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11 // a collaborative effort of two U.S. Department of Energy organizations (Office
12 // of Science and the National Nuclear Security Administration) responsible for
13 // the planning and preparation of a capable exascale ecosystem, including
14 // software, applications, hardware, advanced system engineering and early
15 // testbed platforms, in support of the nation's exascale computing imperative.
16 
17 #include <ceed.h>
18 #include <ceed-backend.h>
19 #include <stdbool.h>
20 #include <stdint.h>
21 #include <string.h>
22 #include "ceed-opt.h"
23 
24 //------------------------------------------------------------------------------
25 // Setup Input/Output Fields
26 //------------------------------------------------------------------------------
27 static int CeedOperatorSetupFields_Opt(CeedQFunction qf, CeedOperator op,
28                                        bool inOrOut, const CeedInt blksize,
29                                        CeedElemRestriction *blkrestr,
30                                        CeedVector *fullevecs, CeedVector *evecs,
31                                        CeedVector *qvecs, CeedInt starte,
32                                        CeedInt numfields, CeedInt Q) {
33   CeedInt dim, ierr, ncomp, size, P;
34   Ceed ceed;
35   ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr);
36   CeedBasis basis;
37   CeedElemRestriction r;
38   CeedOperatorField *opfields;
39   CeedQFunctionField *qffields;
40   if (inOrOut) {
41     ierr = CeedOperatorGetFields(op, NULL, &opfields);
42     CeedChk(ierr);
43     ierr = CeedQFunctionGetFields(qf, NULL, &qffields);
44     CeedChk(ierr);
45   } else {
46     ierr = CeedOperatorGetFields(op, &opfields, NULL);
47     CeedChk(ierr);
48     ierr = CeedQFunctionGetFields(qf, &qffields, NULL);
49     CeedChk(ierr);
50   }
51 
52   // Loop over fields
53   for (CeedInt i=0; i<numfields; i++) {
54     CeedEvalMode emode;
55     ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChk(ierr);
56 
57     if (emode != CEED_EVAL_WEIGHT) {
58       ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &r);
59       CeedChk(ierr);
60       Ceed ceed;
61       ierr = CeedElemRestrictionGetCeed(r, &ceed); CeedChk(ierr);
62       CeedInt nelem, elemsize, lsize, compstride;
63       ierr = CeedElemRestrictionGetNumElements(r, &nelem); CeedChk(ierr);
64       ierr = CeedElemRestrictionGetElementSize(r, &elemsize); CeedChk(ierr);
65       ierr = CeedElemRestrictionGetLVectorSize(r, &lsize); CeedChk(ierr);
66       ierr = CeedElemRestrictionGetNumComponents(r, &ncomp); CeedChk(ierr);
67 
68       bool strided;
69       ierr = CeedElemRestrictionIsStrided(r, &strided); CeedChk(ierr);
70       if (strided) {
71         CeedInt strides[3];
72         ierr = CeedElemRestrictionGetStrides(r, &strides); CeedChk(ierr);
73         ierr = CeedElemRestrictionCreateBlockedStrided(ceed, nelem, elemsize,
74                blksize, ncomp, lsize, strides, &blkrestr[i+starte]);
75         CeedChk(ierr);
76       } else {
77         const CeedInt *offsets = NULL;
78         ierr = CeedElemRestrictionGetOffsets(r, CEED_MEM_HOST, &offsets);
79         CeedChk(ierr);
80         ierr = CeedElemRestrictionGetCompStride(r, &compstride); CeedChk(ierr);
81         ierr = CeedElemRestrictionCreateBlocked(ceed, nelem, elemsize,
82                                                 blksize, ncomp, compstride,
83                                                 lsize, CEED_MEM_HOST,
84                                                 CEED_COPY_VALUES, offsets,
85                                                 &blkrestr[i+starte]);
86         CeedChk(ierr);
87         ierr = CeedElemRestrictionRestoreOffsets(r, &offsets); CeedChk(ierr);
88       }
89       ierr = CeedElemRestrictionCreateVector(blkrestr[i+starte], NULL,
90                                              &fullevecs[i+starte]);
91       CeedChk(ierr);
92     }
93 
94     switch(emode) {
95     case CEED_EVAL_NONE:
96       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr);
97       ierr = CeedVectorCreate(ceed, Q*size*blksize, &evecs[i]); CeedChk(ierr);
98       ierr = CeedVectorCreate(ceed, Q*size*blksize, &qvecs[i]); CeedChk(ierr);
99       break;
100     case CEED_EVAL_INTERP:
101       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr);
102       ierr = CeedElemRestrictionGetElementSize(r, &P);
103       CeedChk(ierr);
104       ierr = CeedVectorCreate(ceed, P*size*blksize, &evecs[i]); CeedChk(ierr);
105       ierr = CeedVectorCreate(ceed, Q*size*blksize, &qvecs[i]); CeedChk(ierr);
106       break;
107     case CEED_EVAL_GRAD:
108       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr);
109       ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr);
110       ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
111       ierr = CeedElemRestrictionGetElementSize(r, &P);
112       CeedChk(ierr);
113       ierr = CeedVectorCreate(ceed, P*size/dim*blksize, &evecs[i]); CeedChk(ierr);
114       ierr = CeedVectorCreate(ceed, Q*size*blksize, &qvecs[i]); CeedChk(ierr);
115       break;
116     case CEED_EVAL_WEIGHT: // Only on input fields
117       ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr);
118       ierr = CeedVectorCreate(ceed, Q*blksize, &qvecs[i]); CeedChk(ierr);
119       ierr = CeedBasisApply(basis, blksize, CEED_NOTRANSPOSE,
120                             CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, qvecs[i]);
121       CeedChk(ierr);
122 
123       break;
124     case CEED_EVAL_DIV:
125       break; // Not implemented
126     case CEED_EVAL_CURL:
127       break; // Not implemented
128     }
129   }
130   return 0;
131 }
132 
133 //------------------------------------------------------------------------------
134 // Setup Operator
135 //------------------------------------------------------------------------------
136 static int CeedOperatorSetup_Opt(CeedOperator op) {
137   int ierr;
138   bool setupdone;
139   ierr = CeedOperatorIsSetupDone(op, &setupdone); CeedChk(ierr);
140   if (setupdone) return 0;
141   Ceed ceed;
142   ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr);
143   Ceed_Opt *ceedimpl;
144   ierr = CeedGetData(ceed, &ceedimpl); CeedChk(ierr);
145   const CeedInt blksize = ceedimpl->blksize;
146   CeedOperator_Opt *impl;
147   ierr = CeedOperatorGetData(op, &impl); CeedChk(ierr);
148   CeedQFunction qf;
149   ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr);
150   CeedInt Q, numinputfields, numoutputfields;
151   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr);
152   ierr = CeedQFunctionIsIdentity(qf, &impl->identityqf); CeedChk(ierr);
153   ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
154   CeedChk(ierr);
155   CeedOperatorField *opinputfields, *opoutputfields;
156   ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields);
157   CeedChk(ierr);
158   CeedQFunctionField *qfinputfields, *qfoutputfields;
159   ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields);
160   CeedChk(ierr);
161 
162   // Allocate
163   ierr = CeedCalloc(numinputfields + numoutputfields, &impl->blkrestr);
164   CeedChk(ierr);
165   ierr = CeedCalloc(numinputfields + numoutputfields, &impl->evecs);
166   CeedChk(ierr);
167   ierr = CeedCalloc(numinputfields + numoutputfields, &impl->edata);
168   CeedChk(ierr);
169 
170   ierr = CeedCalloc(16, &impl->inputstate); CeedChk(ierr);
171   ierr = CeedCalloc(16, &impl->evecsin); CeedChk(ierr);
172   ierr = CeedCalloc(16, &impl->evecsout); CeedChk(ierr);
173   ierr = CeedCalloc(16, &impl->qvecsin); CeedChk(ierr);
174   ierr = CeedCalloc(16, &impl->qvecsout); CeedChk(ierr);
175 
176   impl->numein = numinputfields; impl->numeout = numoutputfields;
177 
178   // Set up infield and outfield pointer arrays
179   // Infields
180   ierr = CeedOperatorSetupFields_Opt(qf, op, 0, blksize, impl->blkrestr,
181                                      impl->evecs, impl->evecsin,
182                                      impl->qvecsin, 0,
183                                      numinputfields, Q);
184   CeedChk(ierr);
185   // Outfields
186   ierr = CeedOperatorSetupFields_Opt(qf, op, 1, blksize, impl->blkrestr,
187                                      impl->evecs, impl->evecsout,
188                                      impl->qvecsout, numinputfields,
189                                      numoutputfields, Q);
190   CeedChk(ierr);
191 
192   // Identity QFunctions
193   if (impl->identityqf) {
194     CeedEvalMode inmode, outmode;
195     CeedQFunctionField *infields, *outfields;
196     ierr = CeedQFunctionGetFields(qf, &infields, &outfields); CeedChk(ierr);
197 
198     for (CeedInt i=0; i<numinputfields; i++) {
199       ierr = CeedQFunctionFieldGetEvalMode(infields[i], &inmode);
200       CeedChk(ierr);
201       ierr = CeedQFunctionFieldGetEvalMode(outfields[i], &outmode);
202       CeedChk(ierr);
203 
204       ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChk(ierr);
205       impl->qvecsout[i] = impl->qvecsin[i];
206       ierr = CeedVectorAddReference(impl->qvecsin[i]); CeedChk(ierr);
207     }
208   }
209 
210   ierr = CeedOperatorSetSetupDone(op); CeedChk(ierr);
211 
212   return 0;
213 }
214 
215 //------------------------------------------------------------------------------
216 // Setup Input Fields
217 //------------------------------------------------------------------------------
218 static inline int CeedOperatorSetupInputs_Opt(CeedInt numinputfields,
219     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
220     CeedVector invec, CeedOperator_Opt *impl, CeedRequest *request) {
221   CeedInt ierr;
222   CeedEvalMode emode;
223   CeedVector vec;
224   uint64_t state;
225 
226   for (CeedInt i=0; i<numinputfields; i++) {
227     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
228     CeedChk(ierr);
229     if (emode == CEED_EVAL_WEIGHT) { // Skip
230     } else {
231       // Get input vector
232       ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr);
233       if (vec != CEED_VECTOR_ACTIVE) {
234         // Restrict
235         ierr = CeedVectorGetState(vec, &state); CeedChk(ierr);
236         if (state != impl->inputstate[i]) {
237           ierr = CeedElemRestrictionApply(impl->blkrestr[i], CEED_NOTRANSPOSE,
238                                           vec, impl->evecs[i], request);
239           CeedChk(ierr);
240           impl->inputstate[i] = state;
241         }
242       } else {
243         // Set Qvec for CEED_EVAL_NONE
244         if (emode == CEED_EVAL_NONE) {
245           ierr = CeedVectorGetArray(impl->evecsin[i], CEED_MEM_HOST,
246                                     &impl->edata[i]); CeedChk(ierr);
247           ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_HOST,
248                                     CEED_USE_POINTER,
249                                     impl->edata[i]); CeedChk(ierr);
250           ierr = CeedVectorRestoreArray(impl->evecsin[i],
251                                         &impl->edata[i]); CeedChk(ierr);
252         }
253       }
254       // Get evec
255       ierr = CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_HOST,
256                                     (const CeedScalar **) &impl->edata[i]);
257       CeedChk(ierr);
258     }
259   }
260   return 0;
261 }
262 
263 //------------------------------------------------------------------------------
264 // Input Basis Action
265 //------------------------------------------------------------------------------
266 static inline int CeedOperatorInputBasis_Opt(CeedInt e, CeedInt Q,
267     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
268     CeedInt numinputfields, CeedInt blksize, CeedVector invec, bool skipactive,
269     CeedOperator_Opt *impl, CeedRequest *request) {
270   CeedInt ierr;
271   CeedInt dim, elemsize, size;
272   CeedElemRestriction Erestrict;
273   CeedEvalMode emode;
274   CeedBasis basis;
275   CeedVector vec;
276 
277   for (CeedInt i=0; i<numinputfields; i++) {
278     ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr);
279     // Skip active input
280     if (skipactive) {
281       if (vec == CEED_VECTOR_ACTIVE)
282         continue;
283     }
284 
285     CeedInt activein = 0;
286     // Get elemsize, emode, size
287     ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict);
288     CeedChk(ierr);
289     ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize);
290     CeedChk(ierr);
291     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
292     CeedChk(ierr);
293     ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChk(ierr);
294     // Restrict block active input
295     if (vec == CEED_VECTOR_ACTIVE) {
296       ierr = CeedElemRestrictionApplyBlock(impl->blkrestr[i], e/blksize,
297                                            CEED_NOTRANSPOSE, invec,
298                                            impl->evecsin[i], request);
299       CeedChk(ierr);
300       activein = 1;
301     }
302     // Basis action
303     switch(emode) {
304     case CEED_EVAL_NONE:
305       if (!activein) {
306         ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_HOST,
307                                   CEED_USE_POINTER,
308                                   &impl->edata[i][e*Q*size]); CeedChk(ierr);
309       }
310       break;
311     case CEED_EVAL_INTERP:
312       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
313       CeedChk(ierr);
314       if (!activein) {
315         ierr = CeedVectorSetArray(impl->evecsin[i], CEED_MEM_HOST,
316                                   CEED_USE_POINTER,
317                                   &impl->edata[i][e*elemsize*size]);
318         CeedChk(ierr);
319       }
320       ierr = CeedBasisApply(basis, blksize, CEED_NOTRANSPOSE,
321                             CEED_EVAL_INTERP, impl->evecsin[i],
322                             impl->qvecsin[i]); CeedChk(ierr);
323       break;
324     case CEED_EVAL_GRAD:
325       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
326       CeedChk(ierr);
327       if (!activein) {
328         ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr);
329         ierr = CeedVectorSetArray(impl->evecsin[i], CEED_MEM_HOST,
330                                   CEED_USE_POINTER,
331                                   &impl->edata[i][e*elemsize*size/dim]);
332         CeedChk(ierr);
333       }
334       ierr = CeedBasisApply(basis, blksize, CEED_NOTRANSPOSE,
335                             CEED_EVAL_GRAD, impl->evecsin[i],
336                             impl->qvecsin[i]); CeedChk(ierr);
337       break;
338     case CEED_EVAL_WEIGHT:
339       break;  // No action
340     // LCOV_EXCL_START
341     case CEED_EVAL_DIV:
342     case CEED_EVAL_CURL: {
343       ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis);
344       CeedChk(ierr);
345       Ceed ceed;
346       ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr);
347       return CeedError(ceed, 1, "Ceed evaluation mode not implemented");
348       // LCOV_EXCL_STOP
349     }
350     }
351   }
352   return 0;
353 }
354 
355 //------------------------------------------------------------------------------
356 // Output Basis Action
357 //------------------------------------------------------------------------------
358 static inline int CeedOperatorOutputBasis_Opt(CeedInt e, CeedInt Q,
359     CeedQFunctionField *qfoutputfields, CeedOperatorField *opoutputfields,
360     CeedInt blksize, CeedInt numinputfields, CeedInt numoutputfields,
361     CeedOperator op, CeedVector outvec, CeedOperator_Opt *impl,
362     CeedRequest *request) {
363   CeedInt ierr;
364   CeedElemRestriction Erestrict;
365   CeedEvalMode emode;
366   CeedBasis basis;
367   CeedVector vec;
368 
369   for (CeedInt i=0; i<numoutputfields; i++) {
370     // Get elemsize, emode, size
371     ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict);
372     CeedChk(ierr);
373     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
374     CeedChk(ierr);
375     // Basis action
376     switch(emode) {
377     case CEED_EVAL_NONE:
378       break; // No action
379     case CEED_EVAL_INTERP:
380       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
381       CeedChk(ierr);
382       ierr = CeedBasisApply(basis, blksize, CEED_TRANSPOSE,
383                             CEED_EVAL_INTERP, impl->qvecsout[i],
384                             impl->evecsout[i]); CeedChk(ierr);
385       break;
386     case CEED_EVAL_GRAD:
387       ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis);
388       CeedChk(ierr);
389       ierr = CeedBasisApply(basis, blksize, CEED_TRANSPOSE,
390                             CEED_EVAL_GRAD, impl->qvecsout[i],
391                             impl->evecsout[i]); CeedChk(ierr);
392       break;
393     // LCOV_EXCL_START
394     case CEED_EVAL_WEIGHT: {
395       Ceed ceed;
396       ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr);
397       return CeedError(ceed, 1, "CEED_EVAL_WEIGHT cannot be an output "
398                        "evaluation mode");
399     }
400     case CEED_EVAL_DIV:
401     case CEED_EVAL_CURL: {
402       Ceed ceed;
403       ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr);
404       return CeedError(ceed, 1, "Ceed evaluation mode not implemented");
405       // LCOV_EXCL_STOP
406     }
407     }
408     // Restrict output block
409     // Get output vector
410     ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); CeedChk(ierr);
411     if (vec == CEED_VECTOR_ACTIVE)
412       vec = outvec;
413     // Restrict
414     ierr = CeedElemRestrictionApplyBlock(impl->blkrestr[i+impl->numein],
415                                          e/blksize, CEED_TRANSPOSE,
416                                          impl->evecsout[i], vec, request);
417     CeedChk(ierr);
418   }
419   return 0;
420 }
421 
422 //------------------------------------------------------------------------------
423 // Restore Input Vectors
424 //------------------------------------------------------------------------------
425 static inline int CeedOperatorRestoreInputs_Opt(CeedInt numinputfields,
426     CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields,
427     CeedOperator_Opt *impl) {
428   CeedInt ierr;
429   CeedEvalMode emode;
430 
431   for (CeedInt i=0; i<numinputfields; i++) {
432     ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode);
433     CeedChk(ierr);
434     if (emode == CEED_EVAL_WEIGHT) { // Skip
435     } else {
436       ierr = CeedVectorRestoreArrayRead(impl->evecs[i],
437                                         (const CeedScalar **) &impl->edata[i]);
438       CeedChk(ierr);
439     }
440   }
441   return 0;
442 }
443 
444 //------------------------------------------------------------------------------
445 // Operator Apply
446 //------------------------------------------------------------------------------
447 static int CeedOperatorApplyAdd_Opt(CeedOperator op, CeedVector invec,
448                                     CeedVector outvec, CeedRequest *request) {
449   int ierr;
450   Ceed ceed;
451   ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr);
452   Ceed_Opt *ceedimpl;
453   ierr = CeedGetData(ceed, &ceedimpl); CeedChk(ierr);
454   CeedInt blksize = ceedimpl->blksize;
455   CeedOperator_Opt *impl;
456   ierr = CeedOperatorGetData(op, &impl); CeedChk(ierr);
457   CeedInt Q, numinputfields, numoutputfields, numelements;
458   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChk(ierr);
459   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr);
460   CeedInt nblks = (numelements/blksize) + !!(numelements%blksize);
461   CeedQFunction qf;
462   ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr);
463   ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
464   CeedChk(ierr);
465   CeedOperatorField *opinputfields, *opoutputfields;
466   ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields);
467   CeedChk(ierr);
468   CeedQFunctionField *qfinputfields, *qfoutputfields;
469   ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields);
470   CeedChk(ierr);
471   CeedEvalMode emode;
472 
473   // Setup
474   ierr = CeedOperatorSetup_Opt(op); CeedChk(ierr);
475 
476   // Input Evecs and Restriction
477   ierr = CeedOperatorSetupInputs_Opt(numinputfields, qfinputfields,
478                                      opinputfields, invec, impl, request);
479   CeedChk(ierr);
480 
481   // Output Lvecs, Evecs, and Qvecs
482   for (CeedInt i=0; i<numoutputfields; i++) {
483     // Set Qvec if needed
484     ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode);
485     CeedChk(ierr);
486     if (emode == CEED_EVAL_NONE) {
487       // Set qvec to single block evec
488       ierr = CeedVectorGetArray(impl->evecsout[i], CEED_MEM_HOST,
489                                 &impl->edata[i + numinputfields]);
490       CeedChk(ierr);
491       ierr = CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_HOST,
492                                 CEED_USE_POINTER,
493                                 impl->edata[i + numinputfields]); CeedChk(ierr);
494       ierr = CeedVectorRestoreArray(impl->evecsout[i],
495                                     &impl->edata[i + numinputfields]);
496       CeedChk(ierr);
497     }
498   }
499 
500   // Loop through elements
501   for (CeedInt e=0; e<nblks*blksize; e+=blksize) {
502     // Input basis apply
503     ierr = CeedOperatorInputBasis_Opt(e, Q, qfinputfields, opinputfields,
504                                       numinputfields, blksize, invec, false,
505                                       impl, request); CeedChk(ierr);
506 
507     // Q function
508     if (!impl->identityqf) {
509       ierr = CeedQFunctionApply(qf, Q*blksize, impl->qvecsin, impl->qvecsout);
510       CeedChk(ierr);
511     }
512 
513     // Output basis apply and restrict
514     ierr = CeedOperatorOutputBasis_Opt(e, Q, qfoutputfields, opoutputfields,
515                                        blksize, numinputfields, numoutputfields,
516                                        op, outvec, impl, request);
517     CeedChk(ierr);
518   }
519 
520   // Restore input arrays
521   ierr = CeedOperatorRestoreInputs_Opt(numinputfields, qfinputfields,
522                                        opinputfields, impl);
523   CeedChk(ierr);
524 
525   return 0;
526 }
527 
528 //------------------------------------------------------------------------------
529 // Assemble Linear QFunction
530 //------------------------------------------------------------------------------
531 static int CeedOperatorLinearAssembleQFunction_Opt(CeedOperator op,
532     CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) {
533   int ierr;
534   Ceed ceed;
535   ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr);
536   Ceed_Opt *ceedimpl;
537   ierr = CeedGetData(ceed, &ceedimpl); CeedChk(ierr);
538   const CeedInt blksize = ceedimpl->blksize;
539   CeedOperator_Opt *impl;
540   ierr = CeedOperatorGetData(op, &impl); CeedChk(ierr);
541   CeedInt Q, numinputfields, numoutputfields, numelements, size;
542   ierr = CeedOperatorGetNumElements(op, &numelements); CeedChk(ierr);
543   ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr);
544   CeedInt nblks = (numelements/blksize) + !!(numelements%blksize);
545   CeedQFunction qf;
546   ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr);
547   ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields);
548   CeedChk(ierr);
549   CeedOperatorField *opinputfields, *opoutputfields;
550   ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields);
551   CeedChk(ierr);
552   CeedQFunctionField *qfinputfields, *qfoutputfields;
553   ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields);
554   CeedChk(ierr);
555   CeedVector vec, lvec;
556   CeedInt numactivein = 0, numactiveout = 0;
557   CeedVector *activein = NULL;
558   CeedScalar *a, *tmp;
559 
560   // Setup
561   ierr = CeedOperatorSetup_Opt(op); CeedChk(ierr);
562 
563   // Check for identity
564   if (impl->identityqf)
565     // LCOV_EXCL_START
566     return CeedError(ceed, 1, "Assembling identity qfunctions not supported");
567   // LCOV_EXCL_STOP
568 
569   // Input Evecs and Restriction
570   ierr = CeedOperatorSetupInputs_Opt(numinputfields, qfinputfields,
571                                      opinputfields, NULL, impl, request);
572   CeedChk(ierr);
573 
574   // Count number of active input fields
575   for (CeedInt i=0; i<numinputfields; i++) {
576     // Get input vector
577     ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr);
578     // Check if active input
579     if (vec == CEED_VECTOR_ACTIVE) {
580       ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChk(ierr);
581       ierr = CeedVectorSetValue(impl->qvecsin[i], 0.0); CeedChk(ierr);
582       ierr = CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_HOST, &tmp);
583       CeedChk(ierr);
584       ierr = CeedRealloc(numactivein + size, &activein); CeedChk(ierr);
585       for (CeedInt field=0; field<size; field++) {
586         ierr = CeedVectorCreate(ceed, Q*blksize, &activein[numactivein+field]);
587         CeedChk(ierr);
588         ierr = CeedVectorSetArray(activein[numactivein+field], CEED_MEM_HOST,
589                                   CEED_USE_POINTER, &tmp[field*Q*blksize]);
590         CeedChk(ierr);
591       }
592       numactivein += size;
593       ierr = CeedVectorRestoreArray(impl->qvecsin[i], &tmp); CeedChk(ierr);
594     }
595   }
596 
597   // Count number of active output fields
598   for (CeedInt i=0; i<numoutputfields; i++) {
599     // Get output vector
600     ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); CeedChk(ierr);
601     // Check if active output
602     if (vec == CEED_VECTOR_ACTIVE) {
603       ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); CeedChk(ierr);
604       numactiveout += size;
605     }
606   }
607 
608   // Check sizes
609   if (!numactivein || !numactiveout)
610     // LCOV_EXCL_START
611     return CeedError(ceed, 1, "Cannot assemble QFunction without active inputs "
612                      "and outputs");
613   // LCOV_EXCL_STOP
614 
615   // Setup lvec
616   ierr = CeedVectorCreate(ceed, nblks*blksize*Q*numactivein*numactiveout,
617                           &lvec); CeedChk(ierr);
618   ierr = CeedVectorGetArray(lvec, CEED_MEM_HOST, &a); CeedChk(ierr);
619 
620   // Create output restriction
621   CeedInt strides[3] = {1, Q, numactivein *numactiveout*Q};
622   ierr = CeedElemRestrictionCreateStrided(ceed, numelements, Q,
623                                           numactivein*numactiveout,
624                                           numactivein*numactiveout*numelements*Q,
625                                           strides, rstr); CeedChk(ierr);
626   // Create assembled vector
627   ierr = CeedVectorCreate(ceed, numelements*Q*numactivein*numactiveout,
628                           assembled); CeedChk(ierr);
629 
630   // Loop through elements
631   for (CeedInt e=0; e<nblks*blksize; e+=blksize) {
632     // Input basis apply
633     ierr = CeedOperatorInputBasis_Opt(e, Q, qfinputfields, opinputfields,
634                                       numinputfields, blksize, NULL, true,
635                                       impl, request); CeedChk(ierr);
636 
637     // Assemble QFunction
638     for (CeedInt in=0; in<numactivein; in++) {
639       // Set Inputs
640       ierr = CeedVectorSetValue(activein[in], 1.0); CeedChk(ierr);
641       if (numactivein > 1) {
642         ierr = CeedVectorSetValue(activein[(in+numactivein-1)%numactivein],
643                                   0.0); CeedChk(ierr);
644       }
645       // Set Outputs
646       for (CeedInt out=0; out<numoutputfields; out++) {
647         // Get output vector
648         ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
649         CeedChk(ierr);
650         // Check if active output
651         if (vec == CEED_VECTOR_ACTIVE) {
652           CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_HOST,
653                              CEED_USE_POINTER, a); CeedChk(ierr);
654           ierr = CeedQFunctionFieldGetSize(qfoutputfields[out], &size);
655           CeedChk(ierr);
656           a += size*Q*blksize; // Advance the pointer by the size of the output
657         }
658       }
659       // Apply QFunction
660       ierr = CeedQFunctionApply(qf, Q*blksize, impl->qvecsin, impl->qvecsout);
661       CeedChk(ierr);
662     }
663   }
664 
665   // Un-set output Qvecs to prevent accidental overwrite of Assembled
666   for (CeedInt out=0; out<numoutputfields; out++) {
667     // Get output vector
668     ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec);
669     CeedChk(ierr);
670     // Check if active output
671     if (vec == CEED_VECTOR_ACTIVE) {
672       CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_HOST, CEED_COPY_VALUES,
673                          NULL); CeedChk(ierr);
674     }
675   }
676 
677   // Restore input arrays
678   ierr = CeedOperatorRestoreInputs_Opt(numinputfields, qfinputfields,
679                                        opinputfields, impl);
680   CeedChk(ierr);
681 
682   // Output blocked restriction
683   ierr = CeedVectorRestoreArray(lvec, &a); CeedChk(ierr);
684   ierr = CeedVectorSetValue(*assembled, 0.0); CeedChk(ierr);
685   CeedElemRestriction blkrstr;
686   ierr = CeedElemRestrictionCreateBlockedStrided(ceed, numelements, Q, blksize,
687          numactivein*numactiveout, numactivein*numactiveout*numelements*Q,
688          strides, &blkrstr); CeedChk(ierr);
689   ierr = CeedElemRestrictionApply(blkrstr, CEED_TRANSPOSE, lvec, *assembled,
690                                   request); CeedChk(ierr);
691 
692   // Cleanup
693   for (CeedInt i=0; i<numactivein; i++) {
694     ierr = CeedVectorDestroy(&activein[i]); CeedChk(ierr);
695   }
696   ierr = CeedFree(&activein); CeedChk(ierr);
697   ierr = CeedVectorDestroy(&lvec); CeedChk(ierr);
698   ierr = CeedElemRestrictionDestroy(&blkrstr); CeedChk(ierr);
699 
700   return 0;
701 }
702 
703 //------------------------------------------------------------------------------
704 // Operator Destroy
705 //------------------------------------------------------------------------------
706 static int CeedOperatorDestroy_Opt(CeedOperator op) {
707   int ierr;
708   CeedOperator_Opt *impl;
709   ierr = CeedOperatorGetData(op, &impl); CeedChk(ierr);
710 
711   for (CeedInt i=0; i<impl->numein+impl->numeout; i++) {
712     ierr = CeedElemRestrictionDestroy(&impl->blkrestr[i]); CeedChk(ierr);
713     ierr = CeedVectorDestroy(&impl->evecs[i]); CeedChk(ierr);
714   }
715   ierr = CeedFree(&impl->blkrestr); CeedChk(ierr);
716   ierr = CeedFree(&impl->evecs); CeedChk(ierr);
717   ierr = CeedFree(&impl->edata); CeedChk(ierr);
718   ierr = CeedFree(&impl->inputstate); CeedChk(ierr);
719 
720   for (CeedInt i=0; i<impl->numein; i++) {
721     ierr = CeedVectorDestroy(&impl->evecsin[i]); CeedChk(ierr);
722     ierr = CeedVectorDestroy(&impl->qvecsin[i]); CeedChk(ierr);
723   }
724   ierr = CeedFree(&impl->evecsin); CeedChk(ierr);
725   ierr = CeedFree(&impl->qvecsin); CeedChk(ierr);
726 
727   for (CeedInt i=0; i<impl->numeout; i++) {
728     ierr = CeedVectorDestroy(&impl->evecsout[i]); CeedChk(ierr);
729     ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChk(ierr);
730   }
731   ierr = CeedFree(&impl->evecsout); CeedChk(ierr);
732   ierr = CeedFree(&impl->qvecsout); CeedChk(ierr);
733 
734   ierr = CeedFree(&impl); CeedChk(ierr);
735   return 0;
736 }
737 
738 //------------------------------------------------------------------------------
739 // Operator Create
740 //------------------------------------------------------------------------------
741 int CeedOperatorCreate_Opt(CeedOperator op) {
742   int ierr;
743   Ceed ceed;
744   ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr);
745   Ceed_Opt *ceedimpl;
746   ierr = CeedGetData(ceed, &ceedimpl); CeedChk(ierr);
747   CeedInt blksize = ceedimpl->blksize;
748   CeedOperator_Opt *impl;
749 
750   ierr = CeedCalloc(1, &impl); CeedChk(ierr);
751   ierr = CeedOperatorSetData(op, impl); CeedChk(ierr);
752 
753   if (blksize != 1 && blksize != 8)
754     // LCOV_EXCL_START
755     return CeedError(ceed, 1, "Opt backend cannot use blocksize: %d", blksize);
756   // LCOV_EXCL_STOP
757 
758   ierr = CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction",
759                                 CeedOperatorLinearAssembleQFunction_Opt);
760   CeedChk(ierr);
761   ierr = CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd",
762                                 CeedOperatorApplyAdd_Opt); CeedChk(ierr);
763   ierr = CeedSetBackendFunction(ceed, "Operator", op, "Destroy",
764                                 CeedOperatorDestroy_Opt); CeedChk(ierr);
765   return 0;
766 }
767 //------------------------------------------------------------------------------
768