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-ref.h" 18 19 //------------------------------------------------------------------------------ 20 // Setup Input/Output Fields 21 //------------------------------------------------------------------------------ 22 static int CeedOperatorSetupFields_Ref(CeedQFunction qf, CeedOperator op, 23 bool inOrOut, 24 CeedVector *fullevecs, CeedVector *evecs, 25 CeedVector *qvecs, CeedInt starte, 26 CeedInt numfields, CeedInt Q) { 27 CeedInt dim, ierr, size, P; 28 Ceed ceed; 29 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 30 CeedBasis basis; 31 CeedElemRestriction Erestrict; 32 CeedOperatorField *opfields; 33 CeedQFunctionField *qffields; 34 if (inOrOut) { 35 ierr = CeedOperatorGetFields(op, NULL, &opfields); CeedChk(ierr); 36 ierr = CeedQFunctionGetFields(qf, NULL, &qffields); CeedChk(ierr); 37 } else { 38 ierr = CeedOperatorGetFields(op, &opfields, NULL); CeedChk(ierr); 39 ierr = CeedQFunctionGetFields(qf, &qffields, NULL); CeedChk(ierr); 40 } 41 42 // Loop over fields 43 for (CeedInt i=0; i<numfields; i++) { 44 CeedEvalMode emode; 45 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChk(ierr); 46 47 if (emode != CEED_EVAL_WEIGHT) { 48 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &Erestrict); 49 CeedChk(ierr); 50 ierr = CeedElemRestrictionCreateVector(Erestrict, NULL, 51 &fullevecs[i+starte]); 52 CeedChk(ierr); 53 } 54 55 switch(emode) { 56 case CEED_EVAL_NONE: 57 ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr); 58 ierr = CeedVectorCreate(ceed, Q*size, &qvecs[i]); CeedChk(ierr); 59 break; 60 case CEED_EVAL_INTERP: 61 ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr); 62 ierr = CeedElemRestrictionGetElementSize(Erestrict, &P); 63 CeedChk(ierr); 64 ierr = CeedVectorCreate(ceed, P*size, &evecs[i]); CeedChk(ierr); 65 ierr = CeedVectorCreate(ceed, Q*size, &qvecs[i]); CeedChk(ierr); 66 break; 67 case CEED_EVAL_GRAD: 68 ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr); 69 ierr = CeedQFunctionFieldGetSize(qffields[i], &size); CeedChk(ierr); 70 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 71 ierr = CeedElemRestrictionGetElementSize(Erestrict, &P); 72 CeedChk(ierr); 73 ierr = CeedVectorCreate(ceed, P*size/dim, &evecs[i]); CeedChk(ierr); 74 ierr = CeedVectorCreate(ceed, Q*size, &qvecs[i]); CeedChk(ierr); 75 break; 76 case CEED_EVAL_WEIGHT: // Only on input fields 77 ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr); 78 ierr = CeedVectorCreate(ceed, Q, &qvecs[i]); CeedChk(ierr); 79 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, 80 CEED_VECTOR_NONE, qvecs[i]); CeedChk(ierr); 81 break; 82 case CEED_EVAL_DIV: 83 break; // Not implemented 84 case CEED_EVAL_CURL: 85 break; // Not implemented 86 } 87 } 88 return 0; 89 } 90 91 //------------------------------------------------------------------------------ 92 // Setup Operator 93 //------------------------------------------------------------------------------/* 94 static int CeedOperatorSetup_Ref(CeedOperator op) { 95 int ierr; 96 bool setupdone; 97 ierr = CeedOperatorGetSetupStatus(op, &setupdone); CeedChk(ierr); 98 if (setupdone) return 0; 99 Ceed ceed; 100 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 101 CeedOperator_Ref *impl; 102 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 103 CeedQFunction qf; 104 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 105 CeedInt Q, numinputfields, numoutputfields; 106 ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr); 107 ierr = CeedQFunctionGetIdentityStatus(qf, &impl->identityqf); CeedChk(ierr); 108 ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 109 CeedChk(ierr); 110 CeedOperatorField *opinputfields, *opoutputfields; 111 ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields); 112 CeedChk(ierr); 113 CeedQFunctionField *qfinputfields, *qfoutputfields; 114 ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields); 115 CeedChk(ierr); 116 117 // Allocate 118 ierr = CeedCalloc(numinputfields + numoutputfields, &impl->evecs); 119 CeedChk(ierr); 120 ierr = CeedCalloc(numinputfields + numoutputfields, &impl->edata); 121 CeedChk(ierr); 122 123 ierr = CeedCalloc(16, &impl->inputstate); CeedChk(ierr); 124 ierr = CeedCalloc(16, &impl->evecsin); CeedChk(ierr); 125 ierr = CeedCalloc(16, &impl->evecsout); CeedChk(ierr); 126 ierr = CeedCalloc(16, &impl->qvecsin); CeedChk(ierr); 127 ierr = CeedCalloc(16, &impl->qvecsout); CeedChk(ierr); 128 129 impl->numein = numinputfields; impl->numeout = numoutputfields; 130 131 // Set up infield and outfield evecs and qvecs 132 // Infields 133 ierr = CeedOperatorSetupFields_Ref(qf, op, 0, impl->evecs, 134 impl->evecsin, impl->qvecsin, 0, 135 numinputfields, Q); 136 CeedChk(ierr); 137 // Outfields 138 ierr = CeedOperatorSetupFields_Ref(qf, op, 1, impl->evecs, 139 impl->evecsout, impl->qvecsout, 140 numinputfields, numoutputfields, Q); 141 CeedChk(ierr); 142 143 // Identity QFunctions 144 if (impl->identityqf) { 145 CeedEvalMode inmode, outmode; 146 CeedQFunctionField *infields, *outfields; 147 ierr = CeedQFunctionGetFields(qf, &infields, &outfields); CeedChk(ierr); 148 149 for (CeedInt i=0; i<numinputfields; i++) { 150 ierr = CeedQFunctionFieldGetEvalMode(infields[i], &inmode); 151 CeedChk(ierr); 152 ierr = CeedQFunctionFieldGetEvalMode(outfields[i], &outmode); 153 CeedChk(ierr); 154 155 ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChk(ierr); 156 impl->qvecsout[i] = impl->qvecsin[i]; 157 ierr = CeedVectorAddReference(impl->qvecsin[i]); CeedChk(ierr); 158 } 159 } 160 161 ierr = CeedOperatorSetSetupDone(op); CeedChk(ierr); 162 163 return 0; 164 } 165 166 //------------------------------------------------------------------------------ 167 // Setup Operator Inputs 168 //------------------------------------------------------------------------------ 169 static inline int CeedOperatorSetupInputs_Ref(CeedInt numinputfields, 170 CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 171 CeedVector invec, const bool skipactive, CeedOperator_Ref *impl, 172 CeedRequest *request) { 173 CeedInt ierr; 174 CeedEvalMode emode; 175 CeedVector vec; 176 CeedElemRestriction Erestrict; 177 uint64_t state; 178 179 for (CeedInt i=0; i<numinputfields; i++) { 180 // Get input vector 181 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 182 if (vec == CEED_VECTOR_ACTIVE) { 183 if (skipactive) 184 continue; 185 else 186 vec = invec; 187 } 188 189 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 190 CeedChk(ierr); 191 // Restrict and Evec 192 if (emode == CEED_EVAL_WEIGHT) { // Skip 193 } else { 194 // Restrict 195 ierr = CeedVectorGetState(vec, &state); CeedChk(ierr); 196 // Skip restriction if input is unchanged 197 if (state != impl->inputstate[i] || vec == invec) { 198 ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 199 CeedChk(ierr); 200 ierr = CeedElemRestrictionApply(Erestrict, CEED_NOTRANSPOSE, vec, 201 impl->evecs[i], request); CeedChk(ierr); 202 impl->inputstate[i] = state; 203 } 204 // Get evec 205 ierr = CeedVectorGetArrayRead(impl->evecs[i], CEED_MEM_HOST, 206 (const CeedScalar **) &impl->edata[i]); 207 CeedChk(ierr); 208 } 209 } 210 return 0; 211 } 212 213 //------------------------------------------------------------------------------ 214 // Input Basis Action 215 //------------------------------------------------------------------------------ 216 static inline int CeedOperatorInputBasis_Ref(CeedInt e, CeedInt Q, 217 CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 218 CeedInt numinputfields, const bool skipactive, CeedOperator_Ref *impl) { 219 CeedInt ierr; 220 CeedInt dim, elemsize, size; 221 CeedElemRestriction Erestrict; 222 CeedEvalMode emode; 223 CeedBasis basis; 224 225 for (CeedInt i=0; i<numinputfields; i++) { 226 // Skip active input 227 if (skipactive) { 228 CeedVector vec; 229 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 230 if (vec == CEED_VECTOR_ACTIVE) 231 continue; 232 } 233 // Get elemsize, emode, size 234 ierr = CeedOperatorFieldGetElemRestriction(opinputfields[i], &Erestrict); 235 CeedChk(ierr); 236 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 237 CeedChk(ierr); 238 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 239 CeedChk(ierr); 240 ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChk(ierr); 241 // Basis action 242 switch(emode) { 243 case CEED_EVAL_NONE: 244 ierr = CeedVectorSetArray(impl->qvecsin[i], CEED_MEM_HOST, 245 CEED_USE_POINTER, 246 &impl->edata[i][e*Q*size]); CeedChk(ierr); 247 break; 248 case CEED_EVAL_INTERP: 249 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChk(ierr); 250 ierr = CeedVectorSetArray(impl->evecsin[i], CEED_MEM_HOST, 251 CEED_USE_POINTER, 252 &impl->edata[i][e*elemsize*size]); 253 CeedChk(ierr); 254 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, 255 CEED_EVAL_INTERP, impl->evecsin[i], 256 impl->qvecsin[i]); CeedChk(ierr); 257 break; 258 case CEED_EVAL_GRAD: 259 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); CeedChk(ierr); 260 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 261 ierr = CeedVectorSetArray(impl->evecsin[i], CEED_MEM_HOST, 262 CEED_USE_POINTER, 263 &impl->edata[i][e*elemsize*size/dim]); 264 CeedChk(ierr); 265 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, 266 CEED_EVAL_GRAD, impl->evecsin[i], 267 impl->qvecsin[i]); CeedChk(ierr); 268 break; 269 case CEED_EVAL_WEIGHT: 270 break; // No action 271 // LCOV_EXCL_START 272 case CEED_EVAL_DIV: 273 case CEED_EVAL_CURL: { 274 ierr = CeedOperatorFieldGetBasis(opinputfields[i], &basis); 275 CeedChk(ierr); 276 Ceed ceed; 277 ierr = CeedBasisGetCeed(basis, &ceed); CeedChk(ierr); 278 return CeedError(ceed, 1, "Ceed evaluation mode not implemented"); 279 break; // Not implemented 280 // LCOV_EXCL_STOP 281 } 282 } 283 } 284 return 0; 285 } 286 287 //------------------------------------------------------------------------------ 288 // Output Basis Action 289 //------------------------------------------------------------------------------ 290 static inline int CeedOperatorOutputBasis_Ref(CeedInt e, CeedInt Q, 291 CeedQFunctionField *qfoutputfields, CeedOperatorField *opoutputfields, 292 CeedInt numinputfields, CeedInt numoutputfields, CeedOperator op, 293 CeedOperator_Ref *impl) { 294 CeedInt ierr; 295 CeedInt dim, elemsize, size; 296 CeedElemRestriction Erestrict; 297 CeedEvalMode emode; 298 CeedBasis basis; 299 300 for (CeedInt i=0; i<numoutputfields; i++) { 301 // Get elemsize, emode, size 302 ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 303 CeedChk(ierr); 304 ierr = CeedElemRestrictionGetElementSize(Erestrict, &elemsize); 305 CeedChk(ierr); 306 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 307 CeedChk(ierr); 308 ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); CeedChk(ierr); 309 // Basis action 310 switch(emode) { 311 case CEED_EVAL_NONE: 312 break; // No action 313 case CEED_EVAL_INTERP: 314 ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); 315 CeedChk(ierr); 316 ierr = CeedVectorSetArray(impl->evecsout[i], CEED_MEM_HOST, 317 CEED_USE_POINTER, 318 &impl->edata[i + numinputfields][e*elemsize*size]); 319 CeedChk(ierr); 320 ierr = CeedBasisApply(basis, 1, CEED_TRANSPOSE, 321 CEED_EVAL_INTERP, impl->qvecsout[i], 322 impl->evecsout[i]); CeedChk(ierr); 323 break; 324 case CEED_EVAL_GRAD: 325 ierr = CeedOperatorFieldGetBasis(opoutputfields[i], &basis); 326 CeedChk(ierr); 327 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 328 ierr = CeedVectorSetArray(impl->evecsout[i], CEED_MEM_HOST, 329 CEED_USE_POINTER, 330 &impl->edata[i + numinputfields][e*elemsize*size/dim]); 331 CeedChk(ierr); 332 ierr = CeedBasisApply(basis, 1, CEED_TRANSPOSE, 333 CEED_EVAL_GRAD, impl->qvecsout[i], 334 impl->evecsout[i]); CeedChk(ierr); 335 break; 336 // LCOV_EXCL_START 337 case CEED_EVAL_WEIGHT: { 338 Ceed ceed; 339 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 340 return CeedError(ceed, 1, "CEED_EVAL_WEIGHT cannot be an output " 341 "evaluation mode"); 342 break; // Should not occur 343 } 344 case CEED_EVAL_DIV: 345 case CEED_EVAL_CURL: { 346 Ceed ceed; 347 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 348 return CeedError(ceed, 1, "Ceed evaluation mode not implemented"); 349 // LCOV_EXCL_STOP 350 break; // Not implemented 351 } 352 } 353 } 354 return 0; 355 } 356 357 //------------------------------------------------------------------------------ 358 // Restore Input Vectors 359 //------------------------------------------------------------------------------ 360 static inline int CeedOperatorRestoreInputs_Ref(CeedInt numinputfields, 361 CeedQFunctionField *qfinputfields, CeedOperatorField *opinputfields, 362 const bool skipactive, CeedOperator_Ref *impl) { 363 CeedInt ierr; 364 CeedEvalMode emode; 365 366 for (CeedInt i=0; i<numinputfields; i++) { 367 // Skip active inputs 368 if (skipactive) { 369 CeedVector vec; 370 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 371 if (vec == CEED_VECTOR_ACTIVE) 372 continue; 373 } 374 // Restore input 375 ierr = CeedQFunctionFieldGetEvalMode(qfinputfields[i], &emode); 376 CeedChk(ierr); 377 if (emode == CEED_EVAL_WEIGHT) { // Skip 378 } else { 379 ierr = CeedVectorRestoreArrayRead(impl->evecs[i], 380 (const CeedScalar **) &impl->edata[i]); 381 CeedChk(ierr); 382 } 383 } 384 return 0; 385 } 386 387 //------------------------------------------------------------------------------ 388 // Operator Apply 389 //------------------------------------------------------------------------------ 390 static int CeedOperatorApply_Ref(CeedOperator op, CeedVector invec, 391 CeedVector outvec, CeedRequest *request) { 392 int ierr; 393 CeedOperator_Ref *impl; 394 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 395 CeedQFunction qf; 396 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 397 CeedInt Q, numelements, numinputfields, numoutputfields, size; 398 ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr); 399 ierr = CeedOperatorGetNumElements(op, &numelements); CeedChk(ierr); 400 ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 401 CeedChk(ierr); 402 CeedOperatorField *opinputfields, *opoutputfields; 403 ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields); 404 CeedChk(ierr); 405 CeedQFunctionField *qfinputfields, *qfoutputfields; 406 ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields); 407 CeedChk(ierr); 408 CeedEvalMode emode; 409 CeedVector vec; 410 CeedElemRestriction Erestrict; 411 412 // Setup 413 ierr = CeedOperatorSetup_Ref(op); CeedChk(ierr); 414 415 // Input Evecs and Restriction 416 ierr = CeedOperatorSetupInputs_Ref(numinputfields, qfinputfields, 417 opinputfields, invec, false, impl, 418 request); CeedChk(ierr); 419 420 // Output Evecs 421 for (CeedInt i=0; i<numoutputfields; i++) { 422 ierr = CeedVectorGetArray(impl->evecs[i+impl->numein], CEED_MEM_HOST, 423 &impl->edata[i + numinputfields]); CeedChk(ierr); 424 } 425 426 // Loop through elements 427 for (CeedInt e=0; e<numelements; e++) { 428 // Output pointers 429 for (CeedInt i=0; i<numoutputfields; i++) { 430 ierr = CeedQFunctionFieldGetEvalMode(qfoutputfields[i], &emode); 431 CeedChk(ierr); 432 if (emode == CEED_EVAL_NONE) { 433 ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); 434 CeedChk(ierr); 435 ierr = CeedVectorSetArray(impl->qvecsout[i], CEED_MEM_HOST, 436 CEED_USE_POINTER, 437 &impl->edata[i + numinputfields][e*Q*size]); 438 CeedChk(ierr); 439 } 440 } 441 442 // Input basis apply 443 ierr = CeedOperatorInputBasis_Ref(e, Q, qfinputfields, opinputfields, 444 numinputfields, false, impl); 445 CeedChk(ierr); 446 447 // Q function 448 if (!impl->identityqf) { 449 ierr = CeedQFunctionApply(qf, Q, impl->qvecsin, impl->qvecsout); 450 CeedChk(ierr); 451 } 452 453 // Output basis apply 454 ierr = CeedOperatorOutputBasis_Ref(e, Q, qfoutputfields, opoutputfields, 455 numinputfields, numoutputfields, op, impl); 456 CeedChk(ierr); 457 } 458 459 // Output restriction 460 for (CeedInt i=0; i<numoutputfields; i++) { 461 // Restore evec 462 ierr = CeedVectorRestoreArray(impl->evecs[i+impl->numein], 463 &impl->edata[i + numinputfields]); 464 CeedChk(ierr); 465 // Get output vector 466 ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); CeedChk(ierr); 467 // Active 468 if (vec == CEED_VECTOR_ACTIVE) 469 vec = outvec; 470 // Restrict 471 ierr = CeedOperatorFieldGetElemRestriction(opoutputfields[i], &Erestrict); 472 CeedChk(ierr); 473 ierr = CeedElemRestrictionApply(Erestrict, CEED_TRANSPOSE, 474 impl->evecs[i+impl->numein], vec, request); 475 CeedChk(ierr); 476 } 477 478 // Restore input arrays 479 ierr = CeedOperatorRestoreInputs_Ref(numinputfields, qfinputfields, 480 opinputfields, false, impl); 481 CeedChk(ierr); 482 483 return 0; 484 } 485 486 //------------------------------------------------------------------------------ 487 // Assemble Linear QFunction 488 //------------------------------------------------------------------------------ 489 static int CeedOperatorAssembleLinearQFunction_Ref(CeedOperator op, 490 CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 491 int ierr; 492 CeedOperator_Ref *impl; 493 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 494 CeedQFunction qf; 495 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 496 CeedInt Q, numelements, numinputfields, numoutputfields, size; 497 ierr = CeedOperatorGetNumQuadraturePoints(op, &Q); CeedChk(ierr); 498 ierr = CeedOperatorGetNumElements(op, &numelements); CeedChk(ierr); 499 ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 500 CeedChk(ierr); 501 CeedOperatorField *opinputfields, *opoutputfields; 502 ierr = CeedOperatorGetFields(op, &opinputfields, &opoutputfields); 503 CeedChk(ierr); 504 CeedQFunctionField *qfinputfields, *qfoutputfields; 505 ierr = CeedQFunctionGetFields(qf, &qfinputfields, &qfoutputfields); 506 CeedChk(ierr); 507 CeedVector vec; 508 CeedInt numactivein = 0, numactiveout = 0; 509 CeedVector *activein = NULL; 510 CeedScalar *a, *tmp; 511 Ceed ceed, ceedparent; 512 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 513 ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent); CeedChk(ierr); 514 ceedparent = ceedparent ? ceedparent : ceed; 515 516 // Setup 517 ierr = CeedOperatorSetup_Ref(op); CeedChk(ierr); 518 519 // Check for identity 520 if (impl->identityqf) 521 // LCOV_EXCL_START 522 return CeedError(ceed, 1, "Assembling identity qfunctions not supported"); 523 // LCOV_EXCL_STOP 524 525 // Input Evecs and Restriction 526 ierr = CeedOperatorSetupInputs_Ref(numinputfields, qfinputfields, 527 opinputfields, NULL, true, impl, request); 528 CeedChk(ierr); 529 530 // Count number of active input fields 531 for (CeedInt i=0; i<numinputfields; i++) { 532 // Get input vector 533 ierr = CeedOperatorFieldGetVector(opinputfields[i], &vec); CeedChk(ierr); 534 // Check if active input 535 if (vec == CEED_VECTOR_ACTIVE) { 536 ierr = CeedQFunctionFieldGetSize(qfinputfields[i], &size); CeedChk(ierr); 537 ierr = CeedVectorSetValue(impl->qvecsin[i], 0.0); CeedChk(ierr); 538 ierr = CeedVectorGetArray(impl->qvecsin[i], CEED_MEM_HOST, &tmp); 539 CeedChk(ierr); 540 ierr = CeedRealloc(numactivein + size, &activein); CeedChk(ierr); 541 for (CeedInt field=0; field<size; field++) { 542 ierr = CeedVectorCreate(ceed, Q, &activein[numactivein+field]); 543 CeedChk(ierr); 544 ierr = CeedVectorSetArray(activein[numactivein+field], CEED_MEM_HOST, 545 CEED_USE_POINTER, &tmp[field*Q]); 546 CeedChk(ierr); 547 } 548 numactivein += size; 549 ierr = CeedVectorRestoreArray(impl->qvecsin[i], &tmp); CeedChk(ierr); 550 } 551 } 552 553 // Count number of active output fields 554 for (CeedInt i=0; i<numoutputfields; i++) { 555 // Get output vector 556 ierr = CeedOperatorFieldGetVector(opoutputfields[i], &vec); CeedChk(ierr); 557 // Check if active output 558 if (vec == CEED_VECTOR_ACTIVE) { 559 ierr = CeedQFunctionFieldGetSize(qfoutputfields[i], &size); CeedChk(ierr); 560 numactiveout += size; 561 } 562 } 563 564 // Check sizes 565 if (!numactivein || !numactiveout) 566 // LCOV_EXCL_START 567 return CeedError(ceed, 1, "Cannot assemble QFunction without active inputs " 568 "and outputs"); 569 // LCOV_EXCL_STOP 570 571 // Create output restriction 572 CeedInt strides[3] = {1, Q, numactivein *numactiveout*Q}; 573 ierr = CeedElemRestrictionCreateStrided(ceedparent, numelements, Q, 574 numelements*Q, 575 numactivein*numactiveout, strides, 576 rstr); CeedChk(ierr); 577 // Create assembled vector 578 ierr = CeedVectorCreate(ceedparent, numelements*Q*numactivein*numactiveout, 579 assembled); CeedChk(ierr); 580 ierr = CeedVectorSetValue(*assembled, 0.0); CeedChk(ierr); 581 ierr = CeedVectorGetArray(*assembled, CEED_MEM_HOST, &a); CeedChk(ierr); 582 583 // Loop through elements 584 for (CeedInt e=0; e<numelements; e++) { 585 // Input basis apply 586 ierr = CeedOperatorInputBasis_Ref(e, Q, qfinputfields, opinputfields, 587 numinputfields, true, impl); 588 CeedChk(ierr); 589 590 // Assemble QFunction 591 for (CeedInt in=0; in<numactivein; in++) { 592 // Set Inputs 593 ierr = CeedVectorSetValue(activein[in], 1.0); CeedChk(ierr); 594 if (numactivein > 1) { 595 ierr = CeedVectorSetValue(activein[(in+numactivein-1)%numactivein], 596 0.0); CeedChk(ierr); 597 } 598 // Set Outputs 599 for (CeedInt out=0; out<numoutputfields; out++) { 600 // Get output vector 601 ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec); 602 CeedChk(ierr); 603 // Check if active output 604 if (vec == CEED_VECTOR_ACTIVE) { 605 CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_HOST, 606 CEED_USE_POINTER, a); CeedChk(ierr); 607 ierr = CeedQFunctionFieldGetSize(qfoutputfields[out], &size); 608 CeedChk(ierr); 609 a += size*Q; // Advance the pointer by the size of the output 610 } 611 } 612 // Apply QFunction 613 ierr = CeedQFunctionApply(qf, Q, impl->qvecsin, impl->qvecsout); 614 CeedChk(ierr); 615 } 616 } 617 618 // Un-set output Qvecs to prevent accidental overwrite of Assembled 619 for (CeedInt out=0; out<numoutputfields; out++) { 620 // Get output vector 621 ierr = CeedOperatorFieldGetVector(opoutputfields[out], &vec); 622 CeedChk(ierr); 623 // Check if active output 624 if (vec == CEED_VECTOR_ACTIVE) { 625 CeedVectorSetArray(impl->qvecsout[out], CEED_MEM_HOST, CEED_COPY_VALUES, 626 NULL); CeedChk(ierr); 627 } 628 } 629 630 // Restore input arrays 631 ierr = CeedOperatorRestoreInputs_Ref(numinputfields, qfinputfields, 632 opinputfields, true, impl); 633 CeedChk(ierr); 634 635 // Restore output 636 ierr = CeedVectorRestoreArray(*assembled, &a); CeedChk(ierr); 637 638 // Cleanup 639 for (CeedInt i=0; i<numactivein; i++) { 640 ierr = CeedVectorDestroy(&activein[i]); CeedChk(ierr); 641 } 642 ierr = CeedFree(&activein); CeedChk(ierr); 643 644 return 0; 645 } 646 647 //------------------------------------------------------------------------------ 648 // Get Basis Emode Pointer 649 //------------------------------------------------------------------------------ 650 static inline void CeedOperatorGetBasisPointer_Ref(CeedScalar **basisptr, 651 CeedEvalMode emode, CeedScalar *identity, CeedScalar *interp, 652 CeedScalar *grad) { 653 switch (emode) { 654 case CEED_EVAL_NONE: 655 *basisptr = identity; 656 break; 657 case CEED_EVAL_INTERP: 658 *basisptr = interp; 659 break; 660 case CEED_EVAL_GRAD: 661 *basisptr = grad; 662 break; 663 case CEED_EVAL_WEIGHT: 664 case CEED_EVAL_DIV: 665 case CEED_EVAL_CURL: 666 break; // Caught by QF Assembly 667 } 668 } 669 670 //------------------------------------------------------------------------------ 671 // Assemble Linear Diagonal 672 //------------------------------------------------------------------------------ 673 static int CeedOperatorAssembleLinearDiagonal_Ref(CeedOperator op, 674 CeedVector *assembled, CeedRequest *request) { 675 int ierr; 676 677 // Assemble QFunction 678 CeedQFunction qf; 679 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 680 CeedInt numinputfields, numoutputfields; 681 ierr= CeedQFunctionGetNumArgs(qf, &numinputfields, &numoutputfields); 682 CeedChk(ierr); 683 CeedVector assembledqf; 684 CeedElemRestriction rstr; 685 ierr = CeedOperatorAssembleLinearQFunction(op, &assembledqf, &rstr, request); 686 CeedChk(ierr); 687 ierr = CeedElemRestrictionDestroy(&rstr); CeedChk(ierr); 688 CeedScalar maxnorm = 0; 689 ierr = CeedVectorNorm(assembledqf, CEED_NORM_MAX, &maxnorm); CeedChk(ierr); 690 691 // Determine active input basis 692 CeedOperatorField *opfields; 693 CeedQFunctionField *qffields; 694 ierr = CeedOperatorGetFields(op, &opfields, NULL); CeedChk(ierr); 695 ierr = CeedQFunctionGetFields(qf, &qffields, NULL); CeedChk(ierr); 696 CeedInt numemodein = 0, ncomp, dim = 1; 697 CeedEvalMode *emodein = NULL; 698 CeedBasis basisin = NULL; 699 CeedElemRestriction rstrin = NULL; 700 for (CeedInt i=0; i<numinputfields; i++) { 701 CeedVector vec; 702 ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChk(ierr); 703 if (vec == CEED_VECTOR_ACTIVE) { 704 ierr = CeedOperatorFieldGetBasis(opfields[i], &basisin); CeedChk(ierr); 705 ierr = CeedBasisGetNumComponents(basisin, &ncomp); CeedChk(ierr); 706 ierr = CeedBasisGetDimension(basisin, &dim); CeedChk(ierr); 707 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstrin); 708 CeedChk(ierr); 709 CeedEvalMode emode; 710 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); 711 CeedChk(ierr); 712 switch (emode) { 713 case CEED_EVAL_NONE: 714 case CEED_EVAL_INTERP: 715 ierr = CeedRealloc(numemodein + 1, &emodein); CeedChk(ierr); 716 emodein[numemodein] = emode; 717 numemodein += 1; 718 break; 719 case CEED_EVAL_GRAD: 720 ierr = CeedRealloc(numemodein + dim, &emodein); CeedChk(ierr); 721 for (CeedInt d=0; d<dim; d++) 722 emodein[numemodein+d] = emode; 723 numemodein += dim; 724 break; 725 case CEED_EVAL_WEIGHT: 726 case CEED_EVAL_DIV: 727 case CEED_EVAL_CURL: 728 break; // Caught by QF Assembly 729 } 730 } 731 } 732 733 // Determine active output basis 734 ierr = CeedOperatorGetFields(op, NULL, &opfields); CeedChk(ierr); 735 ierr = CeedQFunctionGetFields(qf, NULL, &qffields); CeedChk(ierr); 736 CeedInt numemodeout = 0; 737 CeedEvalMode *emodeout = NULL; 738 CeedBasis basisout = NULL; 739 CeedElemRestriction rstrout = NULL; 740 for (CeedInt i=0; i<numoutputfields; i++) { 741 CeedVector vec; 742 ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChk(ierr); 743 if (vec == CEED_VECTOR_ACTIVE) { 744 ierr = CeedOperatorFieldGetBasis(opfields[i], &basisout); CeedChk(ierr); 745 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstrout); 746 CeedChk(ierr); 747 CeedEvalMode emode; 748 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChk(ierr); 749 switch (emode) { 750 case CEED_EVAL_NONE: 751 case CEED_EVAL_INTERP: 752 ierr = CeedRealloc(numemodeout + 1, &emodeout); CeedChk(ierr); 753 emodeout[numemodeout] = emode; 754 numemodeout += 1; 755 break; 756 case CEED_EVAL_GRAD: 757 ierr = CeedRealloc(numemodeout + dim, &emodeout); CeedChk(ierr); 758 for (CeedInt d=0; d<dim; d++) 759 emodeout[numemodeout+d] = emode; 760 numemodeout += dim; 761 break; 762 case CEED_EVAL_WEIGHT: 763 case CEED_EVAL_DIV: 764 case CEED_EVAL_CURL: 765 break; // Caught by QF Assembly 766 } 767 } 768 } 769 770 // Create diagonal vector 771 CeedVector elemdiag; 772 ierr = CeedElemRestrictionCreateVector(rstrin, assembled, &elemdiag); 773 CeedChk(ierr); 774 775 // Assemble element operator diagonals 776 CeedScalar *elemdiagarray, *assembledqfarray; 777 ierr = CeedVectorSetValue(elemdiag, 0.0); CeedChk(ierr); 778 ierr = CeedVectorGetArray(elemdiag, CEED_MEM_HOST, &elemdiagarray); 779 CeedChk(ierr); 780 ierr = CeedVectorGetArray(assembledqf, CEED_MEM_HOST, &assembledqfarray); 781 CeedChk(ierr); 782 CeedInt nelem, nnodes, nqpts; 783 ierr = CeedElemRestrictionGetNumElements(rstrin, &nelem); CeedChk(ierr); 784 ierr = CeedBasisGetNumNodes(basisin, &nnodes); CeedChk(ierr); 785 ierr = CeedBasisGetNumQuadraturePoints(basisin, &nqpts); CeedChk(ierr); 786 // Basis matrices 787 CeedScalar *identity = NULL, *interpin, *interpout, *gradin, *gradout; 788 bool evalNone = false; 789 for (CeedInt i=0; i<numemodein; i++) 790 evalNone = evalNone || (emodein[i] == CEED_EVAL_NONE); 791 for (CeedInt i=0; i<numemodeout; i++) 792 evalNone = evalNone || (emodeout[i] == CEED_EVAL_NONE); 793 if (evalNone) { 794 ierr = CeedCalloc(nqpts*nnodes, &identity); CeedChk(ierr); 795 for (CeedInt i=0; i<(nnodes<nqpts?nnodes:nqpts); i++) 796 identity[i*nnodes+i] = 1.0; 797 } 798 ierr = CeedBasisGetInterp(basisin, &interpin); CeedChk(ierr); 799 ierr = CeedBasisGetInterp(basisout, &interpout); CeedChk(ierr); 800 ierr = CeedBasisGetGrad(basisin, &gradin); CeedChk(ierr); 801 ierr = CeedBasisGetGrad(basisout, &gradout); CeedChk(ierr); 802 // Compute the diagonal of B^T D B 803 // Each element 804 for (CeedInt e=0; e<nelem; e++) { 805 CeedInt dout = -1; 806 // Each basis eval mode pair 807 for (CeedInt eout=0; eout<numemodeout; eout++) { 808 CeedScalar *bt = NULL; 809 if (emodeout[eout] == CEED_EVAL_GRAD) 810 dout += 1; 811 CeedOperatorGetBasisPointer_Ref(&bt, emodeout[eout], identity, interpout, 812 &gradout[dout*nqpts*nnodes]); 813 CeedInt din = -1; 814 for (CeedInt ein=0; ein<numemodein; ein++) { 815 CeedScalar *b = NULL; 816 if (emodein[ein] == CEED_EVAL_GRAD) 817 din += 1; 818 CeedOperatorGetBasisPointer_Ref(&b, emodein[ein], identity, interpin, 819 &gradin[din*nqpts*nnodes]); 820 // Each component 821 for (CeedInt comp=0; comp<ncomp; comp++) 822 // Each qpoint/node pair 823 for (CeedInt q=0; q<nqpts; q++) { 824 const CeedScalar qfvalue = 825 assembledqfarray[((((e*numemodein+ein)*ncomp+comp)* 826 numemodeout+eout)*ncomp+comp)*nqpts+q]; 827 if (fabs(qfvalue) > maxnorm*1e-12) 828 for (CeedInt n=0; n<nnodes; n++) 829 elemdiagarray[(e*ncomp+comp)*nnodes+n] += bt[q*nnodes+n] * 830 qfvalue * b[q*nnodes+n]; 831 } 832 } 833 } 834 } 835 ierr = CeedVectorRestoreArray(elemdiag, &elemdiagarray); CeedChk(ierr); 836 ierr = CeedVectorRestoreArray(assembledqf, &assembledqfarray); CeedChk(ierr); 837 838 // Assemble local operator diagonal 839 ierr = CeedVectorSetValue(*assembled, 0.0); CeedChk(ierr); 840 ierr = CeedElemRestrictionApply(rstrout, CEED_TRANSPOSE, elemdiag, 841 *assembled, request); CeedChk(ierr); 842 843 // Cleanup 844 ierr = CeedVectorDestroy(&assembledqf); CeedChk(ierr); 845 ierr = CeedVectorDestroy(&elemdiag); CeedChk(ierr); 846 ierr = CeedFree(&emodein); CeedChk(ierr); 847 ierr = CeedFree(&emodeout); CeedChk(ierr); 848 ierr = CeedFree(&identity); CeedChk(ierr); 849 850 return 0; 851 } 852 853 //------------------------------------------------------------------------------ 854 // Create FDM Element Inverse 855 //------------------------------------------------------------------------------ 856 int CeedOperatorCreateFDMElementInverse_Ref(CeedOperator op, 857 CeedOperator *fdminv, 858 CeedRequest *request) { 859 int ierr; 860 Ceed ceed, ceedparent; 861 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 862 ierr = CeedGetOperatorFallbackParentCeed(ceed, &ceedparent); CeedChk(ierr); 863 ceedparent = ceedparent ? ceedparent : ceed; 864 CeedQFunction qf; 865 ierr = CeedOperatorGetQFunction(op, &qf); CeedChk(ierr); 866 867 // Determine active input basis 868 bool interp = false, grad = false; 869 CeedBasis basis = NULL; 870 CeedElemRestriction rstr = NULL; 871 CeedOperatorField *opfields; 872 CeedQFunctionField *qffields; 873 ierr = CeedOperatorGetFields(op, &opfields, NULL); CeedChk(ierr); 874 ierr = CeedQFunctionGetFields(qf, &qffields, NULL); CeedChk(ierr); 875 CeedInt numinputfields; 876 ierr = CeedQFunctionGetNumArgs(qf, &numinputfields, NULL); CeedChk(ierr); 877 for (CeedInt i=0; i<numinputfields; i++) { 878 CeedVector vec; 879 ierr = CeedOperatorFieldGetVector(opfields[i], &vec); CeedChk(ierr); 880 if (vec == CEED_VECTOR_ACTIVE) { 881 CeedEvalMode emode; 882 ierr = CeedQFunctionFieldGetEvalMode(qffields[i], &emode); CeedChk(ierr); 883 interp = interp || emode == CEED_EVAL_INTERP; 884 grad = grad || emode == CEED_EVAL_GRAD; 885 ierr = CeedOperatorFieldGetBasis(opfields[i], &basis); CeedChk(ierr); 886 ierr = CeedOperatorFieldGetElemRestriction(opfields[i], &rstr); 887 CeedChk(ierr); 888 } 889 } 890 if (!basis) 891 // LCOV_EXCL_START 892 return CeedError(ceed, 1, "No active field set"); 893 // LCOV_EXCL_STOP 894 CeedInt P1d, Q1d, elemsize, nqpts, dim, ncomp = 1, nelem = 1, nnodes = 1; 895 ierr = CeedBasisGetNumNodes1D(basis, &P1d); CeedChk(ierr); 896 ierr = CeedBasisGetNumNodes(basis, &elemsize); CeedChk(ierr); 897 ierr = CeedBasisGetNumQuadraturePoints1D(basis, &Q1d); CeedChk(ierr); 898 ierr = CeedBasisGetNumQuadraturePoints(basis, &nqpts); CeedChk(ierr); 899 ierr = CeedBasisGetDimension(basis, &dim); CeedChk(ierr); 900 ierr = CeedBasisGetNumComponents(basis, &ncomp); CeedChk(ierr); 901 ierr = CeedElemRestrictionGetNumElements(rstr, &nelem); CeedChk(ierr); 902 ierr = CeedElemRestrictionGetNumNodes(rstr, &nnodes); CeedChk(ierr); 903 904 // Build and diagonalize 1D Mass and Laplacian 905 bool tensorbasis; 906 ierr = CeedBasisGetTensorStatus(basis, &tensorbasis); CeedChk(ierr); 907 if (!tensorbasis) 908 // LCOV_EXCL_START 909 return CeedError(ceed, 1, "FDMElementInverse only supported for tensor " 910 "bases"); 911 // LCOV_EXCL_STOP 912 CeedScalar *work, *mass, *laplace, *x, *x2, *lambda; 913 ierr = CeedMalloc(Q1d*P1d, &work); CeedChk(ierr); 914 ierr = CeedMalloc(P1d*P1d, &mass); CeedChk(ierr); 915 ierr = CeedMalloc(P1d*P1d, &laplace); CeedChk(ierr); 916 ierr = CeedMalloc(P1d*P1d, &x); CeedChk(ierr); 917 ierr = CeedMalloc(P1d*P1d, &x2); CeedChk(ierr); 918 ierr = CeedMalloc(P1d, &lambda); CeedChk(ierr); 919 // -- Mass 920 CeedScalar *interp1d, *grad1d, *qweight1d; 921 ierr = CeedBasisGetInterp1D(basis, &interp1d); CeedChk(ierr); 922 ierr = CeedBasisGetGrad1D(basis, &grad1d); CeedChk(ierr); 923 ierr = CeedBasisGetQWeights(basis, &qweight1d); CeedChk(ierr); 924 for (CeedInt i=0; i<Q1d; i++) 925 for (CeedInt j=0; j<P1d; j++) 926 work[i+j*Q1d] = interp1d[i*P1d+j]*qweight1d[i]; 927 ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)work, 928 (const CeedScalar *)interp1d, mass, P1d, P1d, Q1d); 929 CeedChk(ierr); 930 // -- Laplacian 931 for (CeedInt i=0; i<Q1d; i++) 932 for (CeedInt j=0; j<P1d; j++) 933 work[i+j*Q1d] = grad1d[i*P1d+j]*qweight1d[i]; 934 ierr = CeedMatrixMultiply(ceed, (const CeedScalar *)work, 935 (const CeedScalar *)grad1d, laplace, P1d, P1d, Q1d); 936 CeedChk(ierr); 937 // -- Diagonalize 938 ierr = CeedSimultaneousDiagonalization(ceed, laplace, mass, x, lambda, P1d); 939 CeedChk(ierr); 940 ierr = CeedFree(&work); CeedChk(ierr); 941 ierr = CeedFree(&mass); CeedChk(ierr); 942 ierr = CeedFree(&laplace); CeedChk(ierr); 943 for (CeedInt i=0; i<P1d; i++) 944 for (CeedInt j=0; j<P1d; j++) 945 x2[i+j*P1d] = x[j+i*P1d]; 946 ierr = CeedFree(&x); CeedChk(ierr); 947 948 // Assemble QFunction 949 CeedVector assembled; 950 CeedElemRestriction rstr_qf; 951 ierr = CeedOperatorAssembleLinearQFunction(op, &assembled, &rstr_qf, 952 request); CeedChk(ierr); 953 ierr = CeedElemRestrictionDestroy(&rstr_qf); CeedChk(ierr); 954 CeedScalar maxnorm = 0; 955 ierr = CeedVectorNorm(assembled, CEED_NORM_MAX, &maxnorm); CeedChk(ierr); 956 957 // Calculate element averages 958 CeedInt nfields = ((interp?1:0) + (grad?dim:0))*((interp?1:0) + (grad?dim:0)); 959 CeedScalar *elemavg; 960 const CeedScalar *assembledarray, *qweightsarray; 961 CeedVector qweights; 962 ierr = CeedVectorCreate(ceedparent, nqpts, &qweights); CeedChk(ierr); 963 ierr = CeedBasisApply(basis, 1, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, 964 CEED_VECTOR_NONE, qweights); CeedChk(ierr); 965 ierr = CeedVectorGetArrayRead(assembled, CEED_MEM_HOST, &assembledarray); 966 CeedChk(ierr); 967 ierr = CeedVectorGetArrayRead(qweights, CEED_MEM_HOST, &qweightsarray); 968 CeedChk(ierr); 969 ierr = CeedCalloc(nelem, &elemavg); CeedChk(ierr); 970 for (CeedInt e=0; e<nelem; e++) { 971 CeedInt count = 0; 972 for (CeedInt q=0; q<nqpts; q++) 973 for (CeedInt i=0; i<ncomp*ncomp*nfields; i++) 974 if (fabs(assembledarray[e*nelem*nqpts*ncomp*ncomp*nfields + 975 i*nqpts + q]) > maxnorm*1e-12) { 976 elemavg[e] += assembledarray[e*nelem*nqpts*ncomp*ncomp*nfields + 977 i*nqpts + q] / qweightsarray[q]; 978 count++; 979 } 980 if (count) 981 elemavg[e] /= count; 982 } 983 ierr = CeedVectorRestoreArrayRead(assembled, &assembledarray); CeedChk(ierr); 984 ierr = CeedVectorDestroy(&assembled); CeedChk(ierr); 985 ierr = CeedVectorRestoreArrayRead(qweights, &qweightsarray); CeedChk(ierr); 986 ierr = CeedVectorDestroy(&qweights); CeedChk(ierr); 987 988 // Build FDM diagonal 989 CeedVector qdata; 990 CeedScalar *qdataarray; 991 ierr = CeedVectorCreate(ceedparent, nelem*ncomp*nnodes, &qdata); CeedChk(ierr); 992 ierr = CeedVectorSetArray(qdata, CEED_MEM_HOST, CEED_COPY_VALUES, NULL); 993 CeedChk(ierr); 994 ierr = CeedVectorGetArray(qdata, CEED_MEM_HOST, &qdataarray); CeedChk(ierr); 995 for (CeedInt e=0; e<nelem; e++) 996 for (CeedInt c=0; c<ncomp; c++) 997 for (CeedInt n=0; n<nnodes; n++) { 998 if (interp) 999 qdataarray[(e*ncomp+c)*nnodes+n] = 1; 1000 if (grad) 1001 for (CeedInt d=0; d<dim; d++) { 1002 CeedInt i = (n / CeedIntPow(P1d, d)) % P1d; 1003 qdataarray[(e*ncomp+c)*nnodes+n] += lambda[i]; 1004 } 1005 qdataarray[(e*ncomp+c)*nnodes+n] = 1 / (elemavg[e] * 1006 qdataarray[(e*ncomp+c)*nnodes+n]); 1007 } 1008 ierr = CeedFree(&elemavg); CeedChk(ierr); 1009 ierr = CeedVectorRestoreArray(qdata, &qdataarray); CeedChk(ierr); 1010 1011 // Setup FDM operator 1012 // -- Basis 1013 CeedBasis fdm_basis; 1014 CeedScalar *graddummy, *qrefdummy, *qweightdummy; 1015 ierr = CeedCalloc(P1d*P1d, &graddummy); CeedChk(ierr); 1016 ierr = CeedCalloc(P1d, &qrefdummy); CeedChk(ierr); 1017 ierr = CeedCalloc(P1d, &qweightdummy); CeedChk(ierr); 1018 ierr = CeedBasisCreateTensorH1(ceedparent, dim, ncomp, P1d, P1d, x2, 1019 graddummy, qrefdummy, qweightdummy, 1020 &fdm_basis); CeedChk(ierr); 1021 ierr = CeedFree(&graddummy); CeedChk(ierr); 1022 ierr = CeedFree(&qrefdummy); CeedChk(ierr); 1023 ierr = CeedFree(&qweightdummy); CeedChk(ierr); 1024 ierr = CeedFree(&x2); CeedChk(ierr); 1025 ierr = CeedFree(&lambda); CeedChk(ierr); 1026 1027 // -- Restriction 1028 CeedElemRestriction rstr_i; 1029 CeedInt strides[3] = {1, nnodes, nnodes*ncomp}; 1030 ierr = CeedElemRestrictionCreateStrided(ceedparent, nelem, nnodes, 1031 nnodes*nelem, ncomp, strides, &rstr_i); 1032 CeedChk(ierr); 1033 // -- QFunction 1034 CeedQFunction mass_qf; 1035 ierr = CeedQFunctionCreateInteriorByName(ceedparent, "MassApply", &mass_qf); 1036 CeedChk(ierr); 1037 // -- Operator 1038 ierr = CeedOperatorCreate(ceedparent, mass_qf, NULL, NULL, fdminv); 1039 CeedChk(ierr); 1040 CeedOperatorSetField(*fdminv, "u", rstr_i, fdm_basis, CEED_VECTOR_ACTIVE); 1041 CeedChk(ierr); 1042 CeedOperatorSetField(*fdminv, "qdata", rstr_i, CEED_BASIS_COLLOCATED, qdata); 1043 CeedChk(ierr); 1044 CeedOperatorSetField(*fdminv, "v", rstr_i, fdm_basis, CEED_VECTOR_ACTIVE); 1045 CeedChk(ierr); 1046 1047 // Cleanup 1048 ierr = CeedVectorDestroy(&qdata); CeedChk(ierr); 1049 ierr = CeedBasisDestroy(&fdm_basis); CeedChk(ierr); 1050 ierr = CeedElemRestrictionDestroy(&rstr_i); CeedChk(ierr); 1051 ierr = CeedQFunctionDestroy(&mass_qf); CeedChk(ierr); 1052 1053 return 0; 1054 } 1055 1056 //------------------------------------------------------------------------------ 1057 // Operator Destroy 1058 //------------------------------------------------------------------------------ 1059 static int CeedOperatorDestroy_Ref(CeedOperator op) { 1060 int ierr; 1061 CeedOperator_Ref *impl; 1062 ierr = CeedOperatorGetData(op, (void *)&impl); CeedChk(ierr); 1063 1064 for (CeedInt i=0; i<impl->numein+impl->numeout; i++) { 1065 ierr = CeedVectorDestroy(&impl->evecs[i]); CeedChk(ierr); 1066 } 1067 ierr = CeedFree(&impl->evecs); CeedChk(ierr); 1068 ierr = CeedFree(&impl->edata); CeedChk(ierr); 1069 ierr = CeedFree(&impl->inputstate); CeedChk(ierr); 1070 1071 for (CeedInt i=0; i<impl->numein; i++) { 1072 ierr = CeedVectorDestroy(&impl->evecsin[i]); CeedChk(ierr); 1073 ierr = CeedVectorDestroy(&impl->qvecsin[i]); CeedChk(ierr); 1074 } 1075 ierr = CeedFree(&impl->evecsin); CeedChk(ierr); 1076 ierr = CeedFree(&impl->qvecsin); CeedChk(ierr); 1077 1078 for (CeedInt i=0; i<impl->numeout; i++) { 1079 ierr = CeedVectorDestroy(&impl->evecsout[i]); CeedChk(ierr); 1080 ierr = CeedVectorDestroy(&impl->qvecsout[i]); CeedChk(ierr); 1081 } 1082 ierr = CeedFree(&impl->evecsout); CeedChk(ierr); 1083 ierr = CeedFree(&impl->qvecsout); CeedChk(ierr); 1084 1085 ierr = CeedFree(&impl); CeedChk(ierr); 1086 return 0; 1087 } 1088 1089 //------------------------------------------------------------------------------ 1090 // Operator Create 1091 //------------------------------------------------------------------------------ 1092 int CeedOperatorCreate_Ref(CeedOperator op) { 1093 int ierr; 1094 Ceed ceed; 1095 ierr = CeedOperatorGetCeed(op, &ceed); CeedChk(ierr); 1096 CeedOperator_Ref *impl; 1097 1098 ierr = CeedCalloc(1, &impl); CeedChk(ierr); 1099 ierr = CeedOperatorSetData(op, (void *)&impl); CeedChk(ierr); 1100 1101 ierr = CeedSetBackendFunction(ceed, "Operator", op, "AssembleLinearQFunction", 1102 CeedOperatorAssembleLinearQFunction_Ref); 1103 CeedChk(ierr); 1104 ierr = CeedSetBackendFunction(ceed, "Operator", op, "AssembleLinearDiagonal", 1105 CeedOperatorAssembleLinearDiagonal_Ref); 1106 CeedChk(ierr); 1107 ierr = CeedSetBackendFunction(ceed, "Operator", op, "CreateFDMElementInverse", 1108 CeedOperatorCreateFDMElementInverse_Ref); 1109 CeedChk(ierr); 1110 ierr = CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", 1111 CeedOperatorApply_Ref); CeedChk(ierr); 1112 ierr = CeedSetBackendFunction(ceed, "Operator", op, "Destroy", 1113 CeedOperatorDestroy_Ref); CeedChk(ierr); 1114 return 0; 1115 } 1116 //------------------------------------------------------------------------------ 1117