xref: /petsc/src/dm/dt/fe/impls/vector/fevector.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
12dce792eSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
22dce792eSToby Isaac #include <petsc/private/petscimpl.h>
32dce792eSToby Isaac 
42dce792eSToby Isaac typedef struct _n_PetscFE_Vec {
52dce792eSToby Isaac   PetscFE   scalar_fe;
62dce792eSToby Isaac   PetscInt  num_copies;
72dce792eSToby Isaac   PetscBool interleave_basis;
82dce792eSToby Isaac   PetscBool interleave_components;
92dce792eSToby Isaac } PetscFE_Vec;
102dce792eSToby Isaac 
PetscFEDestroy_Vector(PetscFE fe)112dce792eSToby Isaac static PetscErrorCode PetscFEDestroy_Vector(PetscFE fe)
122dce792eSToby Isaac {
132dce792eSToby Isaac   PetscFE_Vec *v;
142dce792eSToby Isaac 
152dce792eSToby Isaac   PetscFunctionBegin;
162dce792eSToby Isaac   v = (PetscFE_Vec *)fe->data;
172dce792eSToby Isaac   PetscCall(PetscFEDestroy(&v->scalar_fe));
182dce792eSToby Isaac   PetscCall(PetscFree(v));
192dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
202dce792eSToby Isaac }
212dce792eSToby Isaac 
PetscFEView_Vector_Ascii(PetscFE fe,PetscViewer v)222dce792eSToby Isaac static PetscErrorCode PetscFEView_Vector_Ascii(PetscFE fe, PetscViewer v)
232dce792eSToby Isaac {
242dce792eSToby Isaac   PetscInt          dim, Nc, scalar_Nc;
252dce792eSToby Isaac   PetscSpace        basis = NULL;
262dce792eSToby Isaac   PetscDualSpace    dual  = NULL;
272dce792eSToby Isaac   PetscQuadrature   quad  = NULL;
282dce792eSToby Isaac   PetscFE_Vec      *vec;
292dce792eSToby Isaac   PetscViewerFormat fmt;
302dce792eSToby Isaac 
312dce792eSToby Isaac   PetscFunctionBegin;
322dce792eSToby Isaac   vec = (PetscFE_Vec *)fe->data;
332dce792eSToby Isaac   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
342dce792eSToby Isaac   PetscCall(PetscFEGetNumComponents(fe, &Nc));
352dce792eSToby Isaac   PetscCall(PetscFEGetNumComponents(vec->scalar_fe, &scalar_Nc));
362dce792eSToby Isaac   PetscCall(PetscFEGetBasisSpace(fe, &basis));
372dce792eSToby Isaac   PetscCall(PetscFEGetDualSpace(fe, &dual));
382dce792eSToby Isaac   PetscCall(PetscFEGetQuadrature(fe, &quad));
392dce792eSToby Isaac   PetscCall(PetscViewerGetFormat(v, &fmt));
402dce792eSToby Isaac   PetscCall(PetscViewerASCIIPushTab(v));
412dce792eSToby Isaac   if (scalar_Nc == 1) {
422dce792eSToby Isaac     PetscCall(PetscViewerASCIIPrintf(v, "Vector Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components\n", dim, Nc));
432dce792eSToby Isaac   } else {
442dce792eSToby Isaac     PetscCall(PetscViewerASCIIPrintf(v, "Vector Finite Element in %" PetscInt_FMT " dimensions with %" PetscInt_FMT " components (%" PetscInt_FMT " copies of finite element with %" PetscInt_FMT " components)\n", dim, Nc, vec->num_copies, scalar_Nc));
452dce792eSToby Isaac   }
462dce792eSToby Isaac   if (fmt == PETSC_VIEWER_ASCII_INFO_DETAIL) {
472dce792eSToby Isaac     PetscCall(PetscViewerASCIIPrintf(v, "Original finite element:\n"));
482dce792eSToby Isaac     PetscCall(PetscViewerASCIIPushTab(v));
492dce792eSToby Isaac     PetscCall(PetscFEView(vec->scalar_fe, v));
502dce792eSToby Isaac     PetscCall(PetscViewerASCIIPopTab(v));
512dce792eSToby Isaac   }
522dce792eSToby Isaac   if (basis) PetscCall(PetscSpaceView(basis, v));
532dce792eSToby Isaac   if (dual) PetscCall(PetscDualSpaceView(dual, v));
542dce792eSToby Isaac   if (quad) PetscCall(PetscQuadratureView(quad, v));
552dce792eSToby Isaac   PetscCall(PetscViewerASCIIPopTab(v));
562dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
572dce792eSToby Isaac }
582dce792eSToby Isaac 
PetscFEView_Vector(PetscFE fe,PetscViewer v)592dce792eSToby Isaac static PetscErrorCode PetscFEView_Vector(PetscFE fe, PetscViewer v)
602dce792eSToby Isaac {
61*9f196a02SMartin Diehl   PetscBool isascii;
622dce792eSToby Isaac 
632dce792eSToby Isaac   PetscFunctionBegin;
64*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
65*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscFEView_Vector_Ascii(fe, v));
662dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
672dce792eSToby Isaac }
682dce792eSToby Isaac 
PetscFESetUp_Vector(PetscFE fe)692dce792eSToby Isaac static PetscErrorCode PetscFESetUp_Vector(PetscFE fe)
702dce792eSToby Isaac {
712dce792eSToby Isaac   PetscFE_Vec        *v = (PetscFE_Vec *)fe->data;
722dce792eSToby Isaac   PetscDualSpace      dsp;
732dce792eSToby Isaac   PetscInt            n, Ncopies = v->num_copies;
742dce792eSToby Isaac   PetscInt            scalar_n;
752dce792eSToby Isaac   PetscInt           *d, *d_mapped;
762dce792eSToby Isaac   PetscDualSpace_Sum *sum;
772dce792eSToby Isaac   PetscBool           is_sum;
782dce792eSToby Isaac 
792dce792eSToby Isaac   PetscFunctionBegin;
802dce792eSToby Isaac   PetscCall(PetscFESetUp(v->scalar_fe));
812dce792eSToby Isaac   PetscCall(PetscFEGetDimension(v->scalar_fe, &scalar_n));
822dce792eSToby Isaac   PetscCall(PetscFEGetDualSpace(fe, &dsp));
832dce792eSToby Isaac   PetscCall(PetscObjectTypeCompare((PetscObject)dsp, PETSCDUALSPACESUM, &is_sum));
842dce792eSToby Isaac   PetscCheck(is_sum, PetscObjectComm((PetscObject)fe), PETSC_ERR_ARG_INCOMP, "Expected PETSCDUALSPACESUM dual space");
852dce792eSToby Isaac   sum = (PetscDualSpace_Sum *)dsp->data;
862dce792eSToby Isaac   n   = Ncopies * scalar_n;
872dce792eSToby Isaac   PetscCall(PetscCalloc1(n * n, &fe->invV));
882dce792eSToby Isaac   PetscCall(PetscMalloc2(scalar_n, &d, scalar_n, &d_mapped));
892dce792eSToby Isaac   for (PetscInt i = 0; i < scalar_n; i++) d[i] = i;
902dce792eSToby Isaac   for (PetscInt c = 0; c < Ncopies; c++) {
912dce792eSToby Isaac     PetscCall(ISLocalToGlobalMappingApply(sum->all_rows[c], scalar_n, d, d_mapped));
922dce792eSToby Isaac     for (PetscInt i = 0; i < scalar_n; i++) {
932dce792eSToby Isaac       PetscInt         iw      = d_mapped[i];
942dce792eSToby Isaac       PetscReal       *row_w   = &fe->invV[iw * n];
952dce792eSToby Isaac       const PetscReal *row_r   = &v->scalar_fe->invV[i * scalar_n];
962dce792eSToby Isaac       PetscInt         j0      = v->interleave_basis ? c : c * scalar_n;
972dce792eSToby Isaac       PetscInt         jstride = v->interleave_basis ? Ncopies : 1;
982dce792eSToby Isaac 
992dce792eSToby Isaac       for (PetscInt j = 0; j < scalar_n; j++) row_w[j0 + j * jstride] = row_r[j];
1002dce792eSToby Isaac     }
1012dce792eSToby Isaac   }
1022dce792eSToby Isaac   PetscCall(PetscFree2(d, d_mapped));
1032dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1042dce792eSToby Isaac }
1052dce792eSToby Isaac 
PetscFEGetDimension_Vector(PetscFE fe,PetscInt * dim)1062dce792eSToby Isaac static PetscErrorCode PetscFEGetDimension_Vector(PetscFE fe, PetscInt *dim)
1072dce792eSToby Isaac {
1082dce792eSToby Isaac   PetscFE_Vec *v = (PetscFE_Vec *)fe->data;
1092dce792eSToby Isaac 
1102dce792eSToby Isaac   PetscFunctionBegin;
1112dce792eSToby Isaac   PetscCall(PetscFEGetDimension(v->scalar_fe, dim));
1122dce792eSToby Isaac   *dim *= v->num_copies;
1132dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1142dce792eSToby Isaac }
1152dce792eSToby Isaac 
PetscFEVectorInsertTabulation(PetscFE fe,PetscInt npoints,const PetscReal points[],PetscInt k,PetscInt scalar_Nb,PetscInt scalar_point_stride,const PetscReal scalar_Tk[],PetscReal Tk[])1162dce792eSToby Isaac static PetscErrorCode PetscFEVectorInsertTabulation(PetscFE fe, PetscInt npoints, const PetscReal points[], PetscInt k, PetscInt scalar_Nb, PetscInt scalar_point_stride, const PetscReal scalar_Tk[], PetscReal Tk[])
1172dce792eSToby Isaac {
1182dce792eSToby Isaac   PetscFE_Vec *v = (PetscFE_Vec *)fe->data;
1192dce792eSToby Isaac   PetscInt     scalar_Nc;
1202dce792eSToby Isaac   PetscInt     Nc;
1212dce792eSToby Isaac   PetscInt     cdim;
1222dce792eSToby Isaac   PetscInt     dblock;
1232dce792eSToby Isaac   PetscInt     block;
1242dce792eSToby Isaac   PetscInt     scalar_block;
1252dce792eSToby Isaac 
1262dce792eSToby Isaac   PetscFunctionBegin;
1272dce792eSToby Isaac   PetscCall(PetscFEGetNumComponents(v->scalar_fe, &scalar_Nc));
1282dce792eSToby Isaac   Nc = scalar_Nc * v->num_copies;
1292dce792eSToby Isaac   PetscCall(PetscFEGetSpatialDimension(v->scalar_fe, &cdim));
1302dce792eSToby Isaac   dblock       = PetscPowInt(cdim, k);
1312dce792eSToby Isaac   block        = Nc * dblock;
1322dce792eSToby Isaac   scalar_block = scalar_Nc * dblock;
1332dce792eSToby Isaac   for (PetscInt p = 0; p < npoints; p++) {
1342dce792eSToby Isaac     const PetscReal *scalar_Tp = &scalar_Tk[p * scalar_point_stride];
1352dce792eSToby Isaac     PetscReal       *Tp        = &Tk[p * scalar_point_stride * v->num_copies * v->num_copies];
1362dce792eSToby Isaac 
1372dce792eSToby Isaac     for (PetscInt j = 0; j < v->num_copies; j++) {
1382dce792eSToby Isaac       for (PetscInt scalar_i = 0; scalar_i < scalar_Nb; scalar_i++) {
1392dce792eSToby Isaac         PetscInt         i         = v->interleave_basis ? (scalar_i * v->num_copies + j) : (j * scalar_Nb + scalar_i);
1402dce792eSToby Isaac         const PetscReal *scalar_Ti = &scalar_Tp[scalar_i * scalar_block];
1412dce792eSToby Isaac         PetscReal       *Ti        = &Tp[i * block];
1422dce792eSToby Isaac 
1432dce792eSToby Isaac         for (PetscInt scalar_c = 0; scalar_c < scalar_Nc; scalar_c++) {
1442dce792eSToby Isaac           PetscInt         c         = v->interleave_components ? (scalar_c * v->num_copies + j) : (j * scalar_Nc + scalar_c);
1452dce792eSToby Isaac           const PetscReal *scalar_Tc = &scalar_Ti[scalar_c * dblock];
1462dce792eSToby Isaac           PetscReal       *Tc        = &Ti[c * dblock];
1472dce792eSToby Isaac 
1482dce792eSToby Isaac           PetscCall(PetscArraycpy(Tc, scalar_Tc, dblock));
1492dce792eSToby Isaac         }
1502dce792eSToby Isaac       }
1512dce792eSToby Isaac     }
1522dce792eSToby Isaac   }
1532dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1542dce792eSToby Isaac }
1552dce792eSToby Isaac 
PetscFEComputeTabulation_Vector(PetscFE fe,PetscInt npoints,const PetscReal points[],PetscInt K,PetscTabulation T)156ce78bad3SBarry Smith static PetscErrorCode PetscFEComputeTabulation_Vector(PetscFE fe, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation T)
1572dce792eSToby Isaac {
1582dce792eSToby Isaac   PetscFE_Vec    *v = (PetscFE_Vec *)fe->data;
1592dce792eSToby Isaac   PetscInt        scalar_Nc;
1602dce792eSToby Isaac   PetscInt        scalar_Nb;
1612dce792eSToby Isaac   PetscInt        cdim;
1622dce792eSToby Isaac   PetscTabulation scalar_T;
1632dce792eSToby Isaac 
1642dce792eSToby Isaac   PetscFunctionBegin;
1652dce792eSToby Isaac   PetscAssert(npoints == T->Nr * T->Np, PetscObjectComm((PetscObject)fe), PETSC_ERR_PLIB, "Expected to be able decode PetscFECreateTabulation() from PetscTabulation fields");
1662dce792eSToby Isaac   PetscCall(PetscFECreateTabulation(v->scalar_fe, T->Nr, T->Np, points, K, &scalar_T));
1672dce792eSToby Isaac   PetscCall(PetscFEGetNumComponents(v->scalar_fe, &scalar_Nc));
1682dce792eSToby Isaac   PetscCall(PetscFEGetDimension(v->scalar_fe, &scalar_Nb));
1692dce792eSToby Isaac   PetscCall(PetscFEGetSpatialDimension(v->scalar_fe, &cdim));
1702dce792eSToby Isaac   for (PetscInt k = 0; k <= T->K; k++) {
1712dce792eSToby Isaac     PetscReal       *Tk                  = T->T[k];
1722dce792eSToby Isaac     const PetscReal *scalar_Tk           = scalar_T->T[k];
1732dce792eSToby Isaac     PetscInt         dblock              = PetscPowInt(cdim, k);
1742dce792eSToby Isaac     PetscInt         scalar_block        = scalar_Nc * dblock;
1752dce792eSToby Isaac     PetscInt         scalar_point_stride = scalar_Nb * scalar_block;
1762dce792eSToby Isaac 
1772dce792eSToby Isaac     if (v->interleave_basis) {
1782dce792eSToby Isaac       PetscCall(PetscFEVectorInsertTabulation(fe, npoints, points, k, scalar_Nb, scalar_point_stride, scalar_Tk, Tk));
1792dce792eSToby Isaac     } else {
1802dce792eSToby Isaac       PetscDualSpace scalar_dsp;
1812dce792eSToby Isaac       PetscSection   ref_section;
1822dce792eSToby Isaac       PetscInt       pStart, pEnd;
1832dce792eSToby Isaac 
1842dce792eSToby Isaac       PetscCall(PetscFEGetDualSpace(v->scalar_fe, &scalar_dsp));
1852dce792eSToby Isaac       PetscCall(PetscDualSpaceGetSection(scalar_dsp, &ref_section));
1862dce792eSToby Isaac       PetscCall(PetscSectionGetChart(ref_section, &pStart, &pEnd));
1872dce792eSToby Isaac       for (PetscInt p = pStart; p < pEnd; p++) {
1882dce792eSToby Isaac         PetscInt         dof, off;
1892dce792eSToby Isaac         PetscReal       *Tp;
1902dce792eSToby Isaac         const PetscReal *scalar_Tp;
1912dce792eSToby Isaac 
1922dce792eSToby Isaac         PetscCall(PetscSectionGetDof(ref_section, p, &dof));
1932dce792eSToby Isaac         PetscCall(PetscSectionGetOffset(ref_section, p, &off));
1942dce792eSToby Isaac         scalar_Tp = &scalar_Tk[off * scalar_block];
1952dce792eSToby Isaac         Tp        = &Tk[off * scalar_block * v->num_copies * v->num_copies];
1962dce792eSToby Isaac         PetscCall(PetscFEVectorInsertTabulation(fe, npoints, points, k, dof, scalar_point_stride, scalar_Tp, Tp));
1972dce792eSToby Isaac       }
1982dce792eSToby Isaac     }
1992dce792eSToby Isaac   }
2002dce792eSToby Isaac   PetscCall(PetscTabulationDestroy(&scalar_T));
2012dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2022dce792eSToby Isaac }
2032dce792eSToby Isaac 
PetscFECreatePointTrace_Vector(PetscFE fe,PetscInt refPoint,PetscFE * trFE)2042dce792eSToby Isaac static PetscErrorCode PetscFECreatePointTrace_Vector(PetscFE fe, PetscInt refPoint, PetscFE *trFE)
2052dce792eSToby Isaac {
2062dce792eSToby Isaac   PetscFE_Vec *v = (PetscFE_Vec *)fe->data;
2072dce792eSToby Isaac   PetscFE      scalar_trFE;
2082dce792eSToby Isaac   const char  *name;
2092dce792eSToby Isaac 
2102dce792eSToby Isaac   PetscFunctionBegin;
2112dce792eSToby Isaac   PetscCall(PetscFECreatePointTrace(v->scalar_fe, refPoint, &scalar_trFE));
2122dce792eSToby Isaac   PetscCall(PetscFECreateVector(scalar_trFE, v->num_copies, v->interleave_basis, v->interleave_components, trFE));
2132dce792eSToby Isaac   PetscCall(PetscFEDestroy(&scalar_trFE));
2142dce792eSToby Isaac   PetscCall(PetscObjectGetName((PetscObject)fe, &name));
2152dce792eSToby Isaac   if (name) PetscCall(PetscFESetName(*trFE, name));
2162dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2172dce792eSToby Isaac }
2182dce792eSToby Isaac 
2192dce792eSToby Isaac PETSC_INTERN PetscErrorCode PetscFEIntegrate_Basic(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
2202192575eSBarry Smith PETSC_INTERN PetscErrorCode PetscFEIntegrateBd_Basic(PetscDS, PetscInt, PetscBdPointFn *, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
221989fa639SBrad Aagaard PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridResidual_Basic(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
222e3d591f2SMatthew G. Knepley PETSC_INTERN PetscErrorCode PetscFEIntegrateBdJacobian_Basic(PetscDS, PetscWeakForm, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
223989fa639SBrad Aagaard PETSC_INTERN PetscErrorCode PetscFEIntegrateHybridJacobian_Basic(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
2242dce792eSToby Isaac 
PetscFEInitialize_Vector(PetscFE fe)2252dce792eSToby Isaac static PetscErrorCode PetscFEInitialize_Vector(PetscFE fe)
2262dce792eSToby Isaac {
2272dce792eSToby Isaac   PetscFunctionBegin;
2282dce792eSToby Isaac   fe->ops->setfromoptions          = NULL;
2292dce792eSToby Isaac   fe->ops->setup                   = PetscFESetUp_Vector;
2302dce792eSToby Isaac   fe->ops->view                    = PetscFEView_Vector;
2312dce792eSToby Isaac   fe->ops->destroy                 = PetscFEDestroy_Vector;
2322dce792eSToby Isaac   fe->ops->getdimension            = PetscFEGetDimension_Vector;
2332dce792eSToby Isaac   fe->ops->createpointtrace        = PetscFECreatePointTrace_Vector;
234ce78bad3SBarry Smith   fe->ops->computetabulation       = PetscFEComputeTabulation_Vector;
2352dce792eSToby Isaac   fe->ops->integrate               = PetscFEIntegrate_Basic;
2362dce792eSToby Isaac   fe->ops->integratebd             = PetscFEIntegrateBd_Basic;
2372dce792eSToby Isaac   fe->ops->integrateresidual       = PetscFEIntegrateResidual_Basic;
2382dce792eSToby Isaac   fe->ops->integratebdresidual     = PetscFEIntegrateBdResidual_Basic;
2392dce792eSToby Isaac   fe->ops->integratehybridresidual = PetscFEIntegrateHybridResidual_Basic;
2402dce792eSToby Isaac   fe->ops->integratejacobianaction = NULL;
2412dce792eSToby Isaac   fe->ops->integratejacobian       = PetscFEIntegrateJacobian_Basic;
2422dce792eSToby Isaac   fe->ops->integratebdjacobian     = PetscFEIntegrateBdJacobian_Basic;
2432dce792eSToby Isaac   fe->ops->integratehybridjacobian = PetscFEIntegrateHybridJacobian_Basic;
2442dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2452dce792eSToby Isaac }
2462dce792eSToby Isaac 
2472dce792eSToby Isaac /*MC
2482dce792eSToby Isaac   PETSCFEVECTOR = "vector" - A vector-valued `PetscFE` object that is repeated copies
2492dce792eSToby Isaac   of the same underlying finite element.
2502dce792eSToby Isaac 
2512dce792eSToby Isaac   Level: intermediate
2522dce792eSToby Isaac 
2532dce792eSToby Isaac .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`, `PETSCFEBASIC`, `PetscFECreateVector()`
2542dce792eSToby Isaac M*/
PetscFECreate_Vector(PetscFE fe)2552dce792eSToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreate_Vector(PetscFE fe)
2562dce792eSToby Isaac {
2572dce792eSToby Isaac   PetscFE_Vec *v;
2582dce792eSToby Isaac 
2592dce792eSToby Isaac   PetscFunctionBegin;
2602dce792eSToby Isaac   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
2612dce792eSToby Isaac   PetscCall(PetscNew(&v));
2622dce792eSToby Isaac   fe->data = v;
2632dce792eSToby Isaac 
2642dce792eSToby Isaac   PetscCall(PetscFEInitialize_Vector(fe));
2652dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2662dce792eSToby Isaac }
2672dce792eSToby Isaac 
2682dce792eSToby Isaac /*@
2692dce792eSToby Isaac   PetscFECreateVector - Create a vector-valued `PetscFE` from multiple copies of an underlying
2702dce792eSToby Isaac   `PetscFE`.
2712dce792eSToby Isaac 
2722dce792eSToby Isaac   Collective
2732dce792eSToby Isaac 
2742dce792eSToby Isaac   Input Parameters:
2752dce792eSToby Isaac + scalar_fe             - a `PetscFE` finite element
2762dce792eSToby Isaac . num_copies            - a positive integer
2772dce792eSToby Isaac . interleave_basis      - if `PETSC_TRUE`, the first `num_copies` basis vectors
2782dce792eSToby Isaac                           of the output finite element will be copies of the
2792dce792eSToby Isaac                           first basis vector of `scalar_fe`, and so on for the
2802dce792eSToby Isaac                           other basis vectors; otherwise all of the first-copy
2812dce792eSToby Isaac                           basis vectors will come first, followed by all of the
2822dce792eSToby Isaac                           second-copy, and so on.
2832dce792eSToby Isaac - interleave_components - if `PETSC_TRUE`, the first `num_copies` components
2842dce792eSToby Isaac                           of the output finite element will be copies of the
2852dce792eSToby Isaac                           first component of `scalar_fe`, and so on for the
2862dce792eSToby Isaac                           other components; otherwise all of the first-copy
2872dce792eSToby Isaac                           components will come first, followed by all of the
2882dce792eSToby Isaac                           second-copy, and so on.
2892dce792eSToby Isaac 
2902dce792eSToby Isaac   Output Parameter:
2912dce792eSToby Isaac . vector_fe - a `PetscFE` of type `PETSCFEVECTOR` that represent a discretization space with `num_copies` copies of `scalar_fe`
2922dce792eSToby Isaac 
2932dce792eSToby Isaac   Level: intermediate
2942dce792eSToby Isaac 
2952dce792eSToby Isaac .seealso: `PetscFE`, `PetscFEType`, `PetscFECreate()`, `PetscFESetType()`, `PETSCFEBASIC`, `PETSCFEVECTOR`
2962dce792eSToby Isaac @*/
PetscFECreateVector(PetscFE scalar_fe,PetscInt num_copies,PetscBool interleave_basis,PetscBool interleave_components,PetscFE * vector_fe)2972dce792eSToby Isaac PetscErrorCode PetscFECreateVector(PetscFE scalar_fe, PetscInt num_copies, PetscBool interleave_basis, PetscBool interleave_components, PetscFE *vector_fe)
2982dce792eSToby Isaac {
2992dce792eSToby Isaac   MPI_Comm        comm;
3002dce792eSToby Isaac   PetscFE         fe_vec;
3012dce792eSToby Isaac   PetscFE_Vec    *v;
3022dce792eSToby Isaac   PetscInt        scalar_Nc;
3032dce792eSToby Isaac   PetscQuadrature quad;
3042dce792eSToby Isaac 
3052dce792eSToby Isaac   PetscFunctionBegin;
3062dce792eSToby Isaac   PetscValidHeaderSpecific(scalar_fe, PETSCFE_CLASSID, 1);
3072dce792eSToby Isaac   PetscCall(PetscObjectGetComm((PetscObject)scalar_fe, &comm));
3082dce792eSToby Isaac   PetscCheck(num_copies >= 0, comm, PETSC_ERR_ARG_OUTOFRANGE, "Expected positive number of copies, got %" PetscInt_FMT, num_copies);
3092dce792eSToby Isaac   PetscCall(PetscFECreate(comm, vector_fe));
3102dce792eSToby Isaac   fe_vec = *vector_fe;
3112dce792eSToby Isaac   PetscCall(PetscFESetType(fe_vec, PETSCFEVECTOR));
3122dce792eSToby Isaac   v = (PetscFE_Vec *)fe_vec->data;
3132dce792eSToby Isaac   PetscCall(PetscObjectReference((PetscObject)scalar_fe));
3142dce792eSToby Isaac   v->scalar_fe             = scalar_fe;
3152dce792eSToby Isaac   v->num_copies            = num_copies;
3162dce792eSToby Isaac   v->interleave_basis      = interleave_basis;
3172dce792eSToby Isaac   v->interleave_components = interleave_components;
3182dce792eSToby Isaac   PetscCall(PetscFEGetNumComponents(scalar_fe, &scalar_Nc));
3192dce792eSToby Isaac   PetscCall(PetscFESetNumComponents(fe_vec, scalar_Nc * num_copies));
3202dce792eSToby Isaac   PetscCall(PetscFEGetQuadrature(scalar_fe, &quad));
3212dce792eSToby Isaac   PetscCall(PetscFESetQuadrature(fe_vec, quad));
3222dce792eSToby Isaac   PetscCall(PetscFEGetFaceQuadrature(scalar_fe, &quad));
3232dce792eSToby Isaac   PetscCall(PetscFESetFaceQuadrature(fe_vec, quad));
3242dce792eSToby Isaac   {
3252dce792eSToby Isaac     PetscSpace  scalar_sp;
3262dce792eSToby Isaac     PetscSpace *copies;
3272dce792eSToby Isaac     PetscSpace  sp;
3282dce792eSToby Isaac 
3292dce792eSToby Isaac     PetscCall(PetscFEGetBasisSpace(scalar_fe, &scalar_sp));
3302dce792eSToby Isaac     PetscCall(PetscMalloc1(num_copies, &copies));
3312dce792eSToby Isaac     for (PetscInt i = 0; i < num_copies; i++) {
3322dce792eSToby Isaac       PetscCall(PetscObjectReference((PetscObject)scalar_sp));
3332dce792eSToby Isaac       copies[i] = scalar_sp;
3342dce792eSToby Isaac     }
3352dce792eSToby Isaac     PetscCall(PetscSpaceCreateSum(num_copies, copies, PETSC_TRUE, &sp));
3362dce792eSToby Isaac     PetscCall(PetscSpaceSumSetInterleave(sp, interleave_basis, interleave_components));
3372dce792eSToby Isaac     PetscCall(PetscSpaceSetUp(sp));
3382dce792eSToby Isaac     PetscCall(PetscFESetBasisSpace(fe_vec, sp));
3392dce792eSToby Isaac     PetscCall(PetscSpaceDestroy(&sp));
3402dce792eSToby Isaac     for (PetscInt i = 0; i < num_copies; i++) PetscCall(PetscSpaceDestroy(&copies[i]));
3412dce792eSToby Isaac     PetscCall(PetscFree(copies));
3422dce792eSToby Isaac   }
3432dce792eSToby Isaac   {
3442dce792eSToby Isaac     PetscDualSpace  scalar_sp;
3452dce792eSToby Isaac     PetscDualSpace *copies;
3462dce792eSToby Isaac     PetscDualSpace  sp;
3472dce792eSToby Isaac 
3482dce792eSToby Isaac     PetscCall(PetscFEGetDualSpace(scalar_fe, &scalar_sp));
3492dce792eSToby Isaac     PetscCall(PetscMalloc1(num_copies, &copies));
3502dce792eSToby Isaac     for (PetscInt i = 0; i < num_copies; i++) {
3512dce792eSToby Isaac       PetscCall(PetscObjectReference((PetscObject)scalar_sp));
3522dce792eSToby Isaac       copies[i] = scalar_sp;
3532dce792eSToby Isaac     }
3542dce792eSToby Isaac     PetscCall(PetscDualSpaceCreateSum(num_copies, copies, PETSC_TRUE, &sp));
3552dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetInterleave(sp, interleave_basis, interleave_components));
3562dce792eSToby Isaac     PetscCall(PetscDualSpaceSetUp(sp));
3572dce792eSToby Isaac     PetscCall(PetscFESetDualSpace(fe_vec, sp));
3582dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&sp));
3592dce792eSToby Isaac     for (PetscInt i = 0; i < num_copies; i++) PetscCall(PetscDualSpaceDestroy(&copies[i]));
3602dce792eSToby Isaac     PetscCall(PetscFree(copies));
3612dce792eSToby Isaac   }
3622dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
3632dce792eSToby Isaac }
364