/// <summary>
        /// Uses the given edge rule factory to create quadrature rules for
        /// the boundaries of all cells in <paramref name="mask"/>.
        /// </summary>
        /// <param name="mask">
        /// A <b><see cref="CellMask"/></b> containing all cells to be
        /// integrated over.
        /// </param>
        /// <param name="order">
        /// <see cref="IQuadRuleFactory{T}.GetQuadRuleSet"/>
        /// </param>
        /// <returns>
        /// Quadrature rules for the boundary of all cells in
        /// <paramref name="mask"/> in the volume coordinate system.
        /// </returns>
        public IEnumerable <IChunkRulePair <T> > GetQuadRuleSet(ExecutionMask mask, int order)
        {
            if (!(mask is CellMask))
            {
                throw new ArgumentException("Expecting a cell/volume mask.");
            }
            if (mask.MaskType != MaskType.Geometrical)
            {
                throw new ArgumentException("Expecting a geometrical mask.");
            }
            if (!object.ReferenceEquals(this.context, mask.GridData))
            {
                throw new ArgumentException();
            }

            // helper vars
            int D                = this.context.SpatialDimension;
            var CellToEdge       = this.context.iLogicalCells.Cells2Edges;
            var EdgeToCell       = this.context.iLogicalEdges.CellIndices;
            var FaceIndices      = this.context.iGeomEdges.FaceIndices;
            var TrafoIdx         = this.context.iGeomEdges.Edge2CellTrafoIndex;
            var EdgeToCellTrafos = this.context.iGeomEdges.Edge2CellTrafos;
            int NoOfFaces        = RefElement.NoOfFaces;
            var Scalings         = context.iGeomEdges.Edge2CellTrafos_SqrtGramian;
            int J                = this.context.iLogicalCells.NoOfLocalUpdatedCells;

            // output data structure, temporary
            Dictionary <int, T> cellBndRuleMap = new Dictionary <int, T>();

            // Get edge quad rules
            BitArray edgeBitMask = new BitArray(context.iGeomEdges.Count);

            foreach (Chunk chunk in mask)
            {
                foreach (int cell in chunk.Elements)
                {
                    var LocalCellIndexToEdges_cell = CellToEdge[cell];
                    for (int e = LocalCellIndexToEdges_cell.Length - 1; e >= 0; e--)
                    {
                        int edge = Math.Abs(LocalCellIndexToEdges_cell[e]) - 1;
                        edgeBitMask[edge] = true;
                    }
                    cellBndRuleMap.Add(cell, new T()
                    {
                        NumbersOfNodesPerFace = new int[NoOfFaces]
                    });
                }
            }
            IEnumerable <IChunkRulePair <QuadRule> > edgeRuleMap = edgeRuleFactory.GetQuadRuleSet(
                new EdgeMask(context, edgeBitMask, MaskType.Geometrical), order);


            // build cell boundary rule
            var CellBitMask = mask.GetBitMask();

            foreach (var edgeRule in edgeRuleMap)
            {
                QuadRule rule = edgeRule.Rule;
                Debug.Assert(rule.Nodes.GetLength(1) == Math.Max(D - 1, 1));
                int iEdge0 = edgeRule.Chunk.i0;
                int iEdgeE = edgeRule.Chunk.JE;

                for (int iEdge = iEdge0; iEdge < iEdgeE; iEdge++)
                {
                    for (int kk = 0; kk < 2; kk++)   // loop over in and out
                    {
                        int jCell = EdgeToCell[iEdge, kk];
                        if (jCell < 0 || jCell >= J)
                        {
                            continue;
                        }
                        if (!CellBitMask[jCell])
                        {
                            continue;
                        }

                        int    iTrafo = TrafoIdx[iEdge, kk];
                        var    trafo  = EdgeToCellTrafos[iTrafo];
                        int    iFace  = FaceIndices[iEdge, kk];
                        double scl    = Scalings[iTrafo];

                        NodeSet VolumeNodes = new NodeSet(this.RefElement, rule.NoOfNodes, D);
                        trafo.Transform(rule.Nodes, VolumeNodes);
                        VolumeNodes.LockForever();

                        var cellBndRule = cellBndRuleMap[jCell];
                        AddToCellBoundaryRule(cellBndRule, VolumeNodes, rule.Weights, iFace, scl);
                    }
                }
            }


            // return
            // ======
            {
                int[] Cells = cellBndRuleMap.Keys.ToArray();
                Array.Sort(Cells);

                ChunkRulePair <T>[] R = new ChunkRulePair <T> [Cells.Length];
                for (int i = 0; i < R.Length; i++)
                {
                    int jCell = Cells[i];
                    cellBndRuleMap[jCell].Nodes.LockForever();
                    R[i] = new ChunkRulePair <T>(Chunk.GetSingleElementChunk(jCell), cellBndRuleMap[jCell]);
                }


                return(R);
            }
        }
Exemple #2
0
        /// <summary>
        /// Zips two lists of <see cref="IChunkRulePair{T}"/>s (sorted by
        /// <see cref="Chunk"/>) into a single sorted list of
        /// <see cref="IChunkRulePair{T}"/>s. In this process, the items of
        /// <paramref name="B"/> dominate over the items in
        /// <paramref name="A"/>, see <see cref="Merge"/>.
        /// </summary>
        /// <param name="A">
        /// The inferior quadrature rules (i.e., is dominated by
        /// <paramref name="B"/> in case of conflicts)
        /// </param>
        /// <param name="B">
        /// The dominant quadrature rules.
        /// </param>
        /// <returns>
        /// A combination of the list <paramref name="A"/> and
        /// <paramref name="B"/>.
        /// </returns>
        private static List <IChunkRulePair <TQuadRule> > ZipChunkRulePairs(
            List <IChunkRulePair <TQuadRule> > A, List <IChunkRulePair <TQuadRule> > B)
        {
            if (A.Count == 0)
            {
                return(B);
            }
            else if (B.Count == 0)
            {
                return(A);
            }

            var chunkRulePairs = new List <IChunkRulePair <TQuadRule> >(A.Count + B.Count);

            var eA = A.GetEnumerator();
            var eB = B.GetEnumerator();

            eA.MoveNext();
            IChunkRulePair <TQuadRule> currentA = eA.Current;

            while (eB.MoveNext())
            {
                while (currentA != null && currentA.Chunk.i0 <= eB.Current.Chunk.i0)
                {
                    Chunk chunkBeforeCurrentB = new Chunk()
                    {
                        i0  = currentA.Chunk.i0,
                        Len = eB.Current.Chunk.i0 - currentA.Chunk.i0
                    };

                    if (chunkBeforeCurrentB.Len > 0)
                    {
                        chunkRulePairs.Add(new ChunkRulePair <TQuadRule>(
                                               chunkBeforeCurrentB, currentA.Rule));
                    }

                    Chunk chunkAfterCurrentB = new Chunk()
                    {
                        i0  = eB.Current.Chunk.JE,
                        Len = currentA.Chunk.JE - eB.Current.Chunk.JE
                    };

                    if (chunkAfterCurrentB.Len > 0)
                    {
                        currentA = new ChunkRulePair <TQuadRule>(
                            chunkAfterCurrentB, currentA.Rule);
                    }
                    else
                    {
                        if (eA.MoveNext())
                        {
                            currentA = eA.Current;
                        }
                        else
                        {
                            currentA = null;
                        }
                    }
                }

                chunkRulePairs.Add(eB.Current);
            }

            // Current a has not yet been added
            if (currentA != null)
            {
                chunkRulePairs.Add(currentA);
            }

            // B is finished but some items of A are left
            while (eA.MoveNext())
            {
                chunkRulePairs.Add(currentA);
            }

            return(chunkRulePairs);
        }
        /// <summary>
        /// produces an edge quadrature rule
        /// </summary>
        /// <param name="mask">an edge mask</param>
        /// <param name="order">desired order</param>
        /// <returns></returns>
        public IEnumerable <IChunkRulePair <QuadRule> > GetQuadRuleSet(ExecutionMask mask, int order)
        {
            // init & checks
            // =============

            if (!(mask is EdgeMask))
            {
                throw new ArgumentException("Expecting an edge mask.");
            }
#if DEBUG
            var maskBitMask = mask.GetBitMask();
#endif

            var      Edg2Cel     = this.grd.iGeomEdges.CellIndices;
            var      Edg2Fac     = this.grd.iGeomEdges.FaceIndices;
            int      J           = this.grd.Cells.NoOfLocalUpdatedCells;
            QuadRule DefaultRule = this.RefElement.GetQuadratureRule(order);;

            int myIKrfeEdge = this.grd.Edges.EdgeRefElements.IndexOf(this.RefElement, (a, b) => object.ReferenceEquals(a, b));
            if (myIKrfeEdge < 0)
            {
                throw new ApplicationException("fatal error");
            }


            int[] EdgeIndices = mask.ItemEnum.ToArray();
            int   NoEdg       = EdgeIndices.Length;

            // return value
            ChunkRulePair <QuadRule>[] Ret = new ChunkRulePair <QuadRule> [NoEdg];

            // find cells
            // ==========

            BitArray CellBitMask   = new BitArray(J);
            int[]    Cells         = new int[NoEdg]; // mapping: Edge Index --> Cell Index (both geometrical)
            int[]    Faces         = new int[NoEdg]; // mapping: Edge Index --> Face
            BitArray MaxDomainMask = m_maxDomain.GetBitMask();
            {
                for (int i = 0; i < NoEdg; i++)
                {
                    int iEdge = EdgeIndices[i];
                    if (this.grd.Edges.GetRefElementIndex(iEdge) != myIKrfeEdge)
                    {
                        throw new ArgumentException("illegal edge mask");
                    }

                    if (!(grd.Edges.IsEdgeConformalWithCell1(iEdge) || grd.Edges.IsEdgeConformalWithCell2(iEdge)))
                    {
                        throw new NotSupportedException("For an edge that is not conformal with at least one cell, no edge rule can be created from a cell boundary rule.");
                    }

                    int  jCell0 = Edg2Cel[iEdge, 0];
                    int  jCell1 = Edg2Cel[iEdge, 1];
                    bool conf0  = grd.Edges.IsEdgeConformalWithCell1(iEdge);
                    bool conf1  = grd.Edges.IsEdgeConformalWithCell2(iEdge);


                    // this gives no errors for surface elements in 3D
                    bool Allow0 = MaxDomainMask[jCell0];
                    bool Allow1 = (jCell1 >= 0 && jCell1 < J) ? MaxDomainMask[jCell1] : false;


                    // //this is required for MPI parallel calculations
                    //bool Allow0 = true;// AllowedCells[jCell0];
                    //bool Allow1 = (jCell1 >= 0 && jCell1 < J);// ? AllowedCells[jCell1] : false;

                    if (!Allow0 && !Allow1)
                    {
                        // fallback onto default rule, if allowed

                        //if (this.m_DefaultRuleFallbackAllowed) {
                        //    Cells[i] = -1; // by a negative index, we mark that we take the default rule
                        //    Faces[i] = -1;
                        //    Ret[i] = new ChunkRulePair<QuadRule>(Chunk.GetSingleElementChunk(EdgeIndices[i]), DefaultRule);


                        //} else {
                        throw new ArgumentException("unable to find a cell from which the edge rule can be taken.");
                        //}
                    }
                    else
                    {
                        Debug.Assert(Allow0 || Allow1);

                        if (conf0 && Allow0)
                        {
                            // cell 0 is allowed and conformal:
                            // take this, it won't get better

                            CellBitMask[jCell0] = true;
                            Faces[i]            = Edg2Fac[iEdge, 0];
                            Cells[i]            = jCell0;
                        }
                        else if (conf1 && Allow1)
                        {
                            // cell 1 is allowed and conformal:
                            // take this, it won't get better

                            CellBitMask[jCell1] = true;
                            Faces[i]            = Edg2Fac[iEdge, 1];
                            Cells[i]            = jCell1;
                        }
                        else if (Allow0)
                        {
                            // cell 0 is allowed, but NOT conformal

                            CellBitMask[jCell0] = true;
                            Faces[i]            = -Edg2Fac[iEdge, 0]; // by a negative index, we mark a non-conformal edge
                            Cells[i]            = jCell0;
                        }
                        else if (Allow1)
                        {
                            // cell 1 is allowed, but NOT conformal

                            CellBitMask[jCell1] = true;
                            Faces[i]            = -Edg2Fac[iEdge, 1]; // by a negative index, we mark a non-conformal edge
                            Cells[i]            = jCell1;
                        }
                    }
                }
            }


            // get cell boundary rule
            // ======================
            var CellMask = new CellMask(this.grd, CellBitMask);


            IChunkRulePair <CellBoundaryQuadRule>[] cellBndRule = this.m_cellBndQF.GetQuadRuleSet(CellMask, order).ToArray();
            int[] jCell2PairIdx = new int[J];
            for (int i = 0; i < cellBndRule.Length; i++)
            {
                var chk = cellBndRule[i].Chunk;
                for (int jCell = chk.JE - 1; jCell >= chk.i0; jCell--)
                {
                    jCell2PairIdx[jCell] = i + 555;
                }
            }

            int[] iChunk = new int[NoEdg]; // which chunk for which edge?
            for (int i = 0; i < NoEdg; i++)
            {
                if (Cells[i] >= 0)
                {
                    iChunk[i] = jCell2PairIdx[Cells[i]] - 555;
                }
                else
                {
                    iChunk[i] = int.MinValue;
                }
            }

            // build rule
            // ==========
            {
                for (int i = 0; i < NoEdg; i++)   // loop over edges
                //if (MaxDomainMask[Cells[i]] == false)
                //    Debugger.Break();

                //if (Cells[i] >= 0) {
                {
                    var      CellBndR = cellBndRule[iChunk[i]].Rule;
                    QuadRule qrEdge   = null;

                    if (Faces[i] >= 0)
                    {
                        qrEdge = this.CombineQr(null, CellBndR, Faces[i]);
                    }
                    else
                    {
                        throw new NotSupportedException("currently no support for non-conformal edges.");
                    }

                    qrEdge.Nodes.LockForever();
                    qrEdge.Weights.LockForever();
                    Ret[i] = new ChunkRulePair <QuadRule>(Chunk.GetSingleElementChunk(EdgeIndices[i]), qrEdge);
                    //} else {
                    //    Debug.Assert(Ret[i] != null);
                    //}
                }
            }


            // return
            // ======

#if DEBUG
            for (int i = 0; i < Ret.Length; i++)
            {
                Chunk c = Ret[i].Chunk;
                for (int e = c.i0; e < c.JE; e++)
                {
                    Debug.Assert(maskBitMask[e] == true);
                }
            }
#endif

            return(Ret);

            /*
             *
             * var EdgesThatConcernMe = grd.Edges.Edges4RefElement[m_cellBndQF.RefElement.FaceSimplex];
             * EdgeMask edgMask = (mask as EdgeMask).Intersect(EdgesThatConcernMe);
             *
             * var ret = new Dictionary<Chunk, QuadRule>(mask.NoOfItemsLocally);
             #if DEBUG
             * var EdgesOfInterest = edgMask.GetBitMask();
             * var EdgeTouched = (BitArray)EdgesOfInterest.Clone();
             *
             * edgMask.ToTxtFile("Edges-" + grd.MyRank + ".csv", false);
             #endif
             *
             *
             * WeightInbalance = 0.0;
             *
             * // find all cells that are 'touched' by the edge mask
             * // --------------------------------------------------
             * int J = grd.Cells.NoOfLocalUpdatedCells;
             * BitArray TouchedCells = new BitArray(J, false);
             *
             * var Edg2Cell = grd.Edges.CellIndices;
             *
             * int chunkCnt = 0;
             * foreach (var chnk in mask) {
             *  int EE = chnk.JE;
             *  for (int e = chnk.i0; e < EE; e++) {
             *      if (!(grd.Edges.IsEdgeConformalwithCell1(e) || grd.Edges.IsEdgeConformalwithCell2(e))) {
             *          throw new NotSupportedException("For an edge that is not conformal with at least one cell, no edge rule can be created from a cell boundary rule.");
             *      }
             *
             *      int j1 = Edg2Cell[e, 0], j2 = Edg2Cell[e, 1];
             *      TouchedCells[j1] = true;
             *      if (j2 >= 0 && j2 < J)
             *          TouchedCells[j2] = true;
             *
             *
             *      Chunk singleEdgeChunk;
             *      singleEdgeChunk.i0 = e;
             *      singleEdgeChunk.Len = 1;
             *
             *      ret.Add(singleEdgeChunk, null);
             *  }
             *  chunkCnt++;
             * }
             *
             * CellMask celMask = new CellMask(grd, TouchedCells);
             * celMask.ToTxtFile("CellsU-" + grd.MyRank + ".csv", false);
             * celMask = celMask.Intersect(m_maxDomain);
             * celMask.ToTxtFile("Cells-" + grd.MyRank + ".csv", false);
             *
             * // create cell boundary rule!
             * IEnumerable<IChunkRulePair<CellBoundaryQuadRule>> cellBndRule = m_cellBndQF.GetQuadRuleSet(celMask, order);
             *
             * //// do MPI communication (get rules for external cells)
             * //{
             * //    int size, rank;
             * //    csMPI.Raw.Comm_Rank(csMPI.Raw._COMM.WORLD, out rank);
             * //    csMPI.Raw.Comm_Size(csMPI.Raw._COMM.WORLD, out size);
             *
             * //    if (size > 1)
             * //        throw new NotSupportedException("currently no MPI support");
             * //}
             *
             * // assign the cell boundary rule to edges
             * // --------------------------------------
             * var volSplx = m_cellBndQF.RefElement;
             * int NoOfFaces = volSplx.NoOfFaces;
             * var Cells2Edge = grd.Cells.Cells2Edges;
             * var FaceIndices = grd.Edges.FaceIndices;
             *
             * int cnt = -1;
             * foreach (var kv in cellBndRule) { // loop over cell chunks (in the cell boundary rule)...
             *  Chunk chk = kv.Chunk;
             *  CellBoundaryQuadRule qr = kv.Rule;
             *  cnt++;
             *
             *  int JE = chk.JE;
             *  for (int j = chk.i0; j < JE; j++) { // loop over all cells in chunk...
             *      Debug.Assert(qr.NumbersOfNodesPerFace.Length == NoOfFaces);
             *      var Cells2Edge_j = Cells2Edge[j];
             *
             *      for (int _e = Cells2Edge_j.Length - 1; _e >= 0; _e--) { // loop over all edges of cell...
             *
             *          //if (qr.NumbersOfNodesPerEdge[e] <= 0)
             *          //    // no contribution from this edge
             *          //    continue;
             *
             *          bool isInCell = Cells2Edge_j[_e] >= 0;
             *          int iEdge = Math.Abs(Cells2Edge_j[_e]) - 1;
             *          int iFace = FaceIndices[iEdge, isInCell ? 0 : 1];
             *
             *          if (isInCell) {
             *              if (!grd.Edges.IsEdgeConformalwithCell1(iEdge))
             *                  // we already tested that at least one cell is conformal with the edge,
             *                  // so it is safe to drop this face
             *                  continue;
             *          } else {
             *              if (!grd.Edges.IsEdgeConformalwithCell2(iEdge))
             *                  // we already tested that at least one cell is conformal with the edge,
             *                  // so it is safe to drop this face
             *                  continue;
             *          }
             *
             *
             *          Chunk singleEdgeChunk = Chunk.GetSingleElementChunk(iEdge);
             *
             *          QuadRule qrEdge;
             *          if (ret.TryGetValue(singleEdgeChunk, out qrEdge)) {
             *              // we are interested in this edge!
             *
             #if DEBUG
             *              Debug.Assert(EdgesOfInterest[iEdge] == true);
             *              EdgeTouched[iEdge] = false;
             *
             *              var vtx = this.RefElement.Vertices;
             *              MultidimensionalArray _vtx = MultidimensionalArray.Create(vtx.GetLength(0), vtx.GetLength(1));
             *              _vtx.Set(vtx);
             *
             *              var RefCoord = MultidimensionalArray.Create(vtx.GetLength(0), vtx.GetLength(1) + 1);
             *              var PhysCoord = MultidimensionalArray.Create(1, vtx.GetLength(0), vtx.GetLength(1) + 1);
             *
             *              volSplx.TransformFaceCoordinates(iFace, _vtx, RefCoord);
             *              grd.TransformLocal2Global(RefCoord, PhysCoord, j, 1, 0);
             *
             #endif
             *
             *              qrEdge = CombineQr(qrEdge, qr, iFace, j);
             *
             *              Debug.Assert(qrEdge != null);
             *              ret[singleEdgeChunk] = qrEdge;
             *          } else {
             *              // nop: the edge is not in the 'edgMask'!
             *              continue;
             *          }
             *      }
             *
             *  }
             * }
             *
             #if DEBUG
             * (new EdgeMask(this.grd, EdgeTouched)).ToTxtFile("failedEdges-" + grd.MyRank + ".csv", false);
             *
             * for (int i = EdgeTouched.Length - 1; i >= 0; i--)
             *  Debug.Assert(EdgeTouched[i] == false);
             #endif
             *
             * return ret.Select(p => new ChunkRulePair<QuadRule>(p.Key, p.Value));
             */
        }
Exemple #4
0
        static List <IChunkRulePair <TQuadRule> > ZipHelper2(List <IChunkRulePair <TQuadRule> > A, List <IChunkRulePair <TQuadRule> > B)
        {
            // spezialfaelle
            // -------------
            if (A.Count == 0)
            {
                return(B);
            }
            if (B.Count == 0)
            {
                return(A);
            }

            List <IChunkRulePair <TQuadRule> > output = new List <IChunkRulePair <TQuadRule> >();

            var eA = A.GetEnumerator();
            var eB = B.GetEnumerator();

            bool runA = eA.MoveNext();
            bool runB = eB.MoveNext();
            ChunkRulePair <TQuadRule> cA = null;
            ChunkRulePair <TQuadRule> cB = null;

            if (runA)
            {
                cA = new ChunkRulePair <TQuadRule>(eA.Current);
            }
            if (runB)
            {
                cB = new ChunkRulePair <TQuadRule>(eB.Current);
            }

            // solong in boade lischtn wos drinnen isch,
            // miassn ma se zommentian...
            // -----------------------------------------
            while (runA && runB)
            {
                if (cA._Chunk.JE <= cB._Chunk.i0)
                {
                    // Fall 1
                    output.Add(cA);
                    runA = eA.MoveNext();
                    if (runA)
                    {
                        cA = new ChunkRulePair <TQuadRule>(eA.Current);
                    }
                    continue;
                }

                if (cB._Chunk.JE <= cA._Chunk.i0)
                {
                    output.Add(cB);
                    runB = eB.MoveNext();
                    if (runB)
                    {
                        cB = new ChunkRulePair <TQuadRule>(eB.Current);
                    }
                    continue;
                }


                if (cB.Chunk.i0 <= cA.Chunk.i0)
                {
                    // iatz kimmt zersch't amol cB!
                    output.Add(cB);

                    if (cB._Chunk.JE < cA._Chunk.JE)
                    {
                        // ein stückl von cA bleibt ibrig.
                        // -> aufkalten firn negscht'n loop
                        cA._Chunk.Len = cA._Chunk.JE - cB._Chunk.JE;
                        cA._Chunk.i0  = cB._Chunk.JE;
                        Debug.Assert(cA.Chunk.Len > 0);
                    }
                    else
                    {
                        while (runA && (cB._Chunk.JE >= cA._Chunk.JE))
                        {
                            // cA ist vollständig von cB überdeckt
                            // aktuelles cA wegschmeissen, weiter zum negscht'n
                            runA = eA.MoveNext();
                            if (runA)
                            {
                                cA = new ChunkRulePair <TQuadRule>(eA.Current);
                            }
                        }

                        if (runA && cA._Chunk.i0 <= cB._Chunk.JE)
                        {
                            // vo dem cA wo ma jetz ham, bleibt a stückl übrig.
                            // -> aufkalten fir negscht'n loop
                            cA._Chunk.Len = cA._Chunk.JE - cB._Chunk.JE;
                            cA._Chunk.i0  = cB._Chunk.JE;
                            Debug.Assert(cA.Chunk.Len > 0);
                        }
                    }

                    // weiter zum negscht'n cB
                    runB = eB.MoveNext();
                    if (runB)
                    {
                        cB = new ChunkRulePair <TQuadRule>(eB.Current);
                    }
                    continue;
                }

                if (cB.Chunk.i0 > cA.Chunk.i0)
                {
                    // vorne isch no a stückl cA!
                    Debug.Assert(cA._Chunk.JE > cB._Chunk.i0);
                    Debug.Assert(cB._Chunk.i0 < cA._Chunk.JE);
                    var cX = new ChunkRulePair <TQuadRule>(cA);
                    cX._Chunk.Len = cB._Chunk.i0 - cA._Chunk.i0;
                    Debug.Assert(cX.Chunk.Len > 0);
                    output.Add(cX);

                    cA._Chunk.Len = cA._Chunk.JE - cX._Chunk.JE;
                    Debug.Assert(cA.Chunk.Len > 0); // müsste eigentlich schon durch "Fall 1" (siehe oben) abgefangen werden.
                    cA._Chunk.i0 = cX._Chunk.JE;
                    continue;
                }


                Debug.Assert(false); // wenn mir do her kemmen, nor isch eppes gonz tschelawenget.
            }

            // izant isch nur mehr in A oder in B wos drinnen,
            // oder in koaner vu boaden.
            // -----------------------------------------------

            while (runB)
            {
                output.Add(cB);
                runB = eB.MoveNext();
                if (runB)
                {
                    cB = new ChunkRulePair <TQuadRule>(eB.Current);
                }
            }

            while (runA)
            {
                output.Add(cA);
                runA = eA.MoveNext();
                if (runA)
                {
                    cA = new ChunkRulePair <TQuadRule>(eA.Current);
                }
            }

            // return
            // ------

            return(output);
        }