/// <summary>
        /// Run this
        /// </summary>
        /// <param name="mask"></param>
        /// <param name="order"></param>
        void CalculateComboQuadRuleSet(ExecutionMask mask, int order)
        {
            comboRule.order = order;
            rulez[0].Clear();
            rulez[1].Clear();
            //Find quadrature nodes and weights in each cell/chunk
#if LOG_TIME
            Stopwatch stopWatch = new Stopwatch();
            stopWatch.Start();
#endif
            foreach (Chunk chunk in mask)
            {
                foreach (int cell in chunk.Elements)
                {
                    QuadRule[] sayeRule = comboRule.ComboEvaluate(cell);
                    ChunkRulePair <QuadRule> sayePair_volume  = new ChunkRulePair <QuadRule>(Chunk.GetSingleElementChunk(cell), sayeRule[0]);
                    ChunkRulePair <QuadRule> sayePair_surface = new ChunkRulePair <QuadRule>(Chunk.GetSingleElementChunk(cell), sayeRule[1]);
                    rulez[0].Add(sayePair_volume);
                    rulez[1].Add(sayePair_surface);
                }
            }
#if LOG_TIME
            stopWatch.Stop();
            long ts = stopWatch.ElapsedMilliseconds;
            Console.WriteLine("Calculated combo cutcell rule : {0}ms", ts);
#endif
        }
        public IEnumerable <IChunkRulePair <QuadRule> > GetQuadRuleSet(ExecutionMask mask, int order)
        {
            rule.order = order;
            var result = new List <ChunkRulePair <QuadRule> >();


#if LOG_TIME
            Stopwatch stopWatch = new Stopwatch();
            stopWatch.Start();
#endif
            //Find quadrature nodes and weights in each cell/chunk
            foreach (Chunk chunk in mask)
            {
                foreach (int cell in chunk.Elements)
                {
                    QuadRule sayeRule = rule.Evaluate(cell);
                    ChunkRulePair <QuadRule> sayePair = new ChunkRulePair <QuadRule>(Chunk.GetSingleElementChunk(cell), sayeRule);
                    result.Add(sayePair);
                }
            }
#if LOG_TIME
            stopWatch.Stop();
            long ts = stopWatch.ElapsedMilliseconds;
            Console.WriteLine("Calculated cutcell rule : {0}ms", ts);
#endif
            return(result);
        }
Beispiel #3
0
            public IEnumerable <IChunkRulePair <QuadRule> > GetQuadRuleSet(ExecutionMask mask, int order)
            {
                if (!(mask is CellMask))
                {
                    throw new ArgumentException("Expecting a cell mask.");
                }
                if (mask.MaskType != MaskType.Geometrical)
                {
                    throw new ArgumentException("Expecting a geometrical mask.");
                }

#if DEBUG
                if (mask.Except(m_Owner.MaxGrid).NoOfItemsLocally > 0)
                {
                    throw new NotSupportedException("'mask' must be a subset of the cut cells, for my reference element.");
                }
#endif
                if (!Rules.ContainsKey(order))
                {
                    m_Owner.GetQuadRuleSet_Internal(order);
                }

                if (mask.NoOfItemsLocally == m_Owner.MaxGrid.NoOfItemsLocally)
                {
                    // aggressive
                    return(Rules[order]);
                }
                else
                {
                    var Rule = Rules[order];

                    SubGrid S          = new SubGrid((CellMask)mask);
                    var     jsub2jcell = S.SubgridIndex2LocalCellIndex;
                    var     Ret        = new ChunkRulePair <QuadRule> [S.LocalNoOfCells_WithExternal];

                    int L = Ret.Length, H = Rule.Length;
                    int h = 0;
                    for (int jsub = 0; jsub < L; jsub++)
                    {
                        int jCell = jsub2jcell[jsub];

                        Debug.Assert(Rule[h].Chunk.Len == 1);

                        while (jCell > Rule[h].Chunk.i0)
                        {
                            h++;
                        }

                        Debug.Assert(jCell == Rule[h].Chunk.i0);
                        Ret[jsub] = Rule[h];
                    }

                    return(Ret);
                }
            }
Beispiel #4
0
        void GetQuadRuleSet_Internal(int order)
        {
            using (new FuncTrace()) {
                int original_order = order;
                stpwGetQuadRuleSet.Start();
                order = Math.Max(order, 2); // Ansatz won't work below 2!
                int D = this.tracker.GridDat.SpatialDimension;;

                // check arguments, init
                // =====================

                CellMask _mask = this.MaxGrid;

                // subgrid on which the volume rule should be constructed
                // ======================================================
                CellBoundaryQuadratureScheme cellBndSchme     = new CellBoundaryQuadratureScheme(this.cellFaceFactory, _mask);
                CellBoundaryQuadratureScheme cellBndLineSchme = null;
                if (this.UseStokes)
                {
                    cellBndLineSchme = new CellBoundaryQuadratureScheme(this.LevelSetBoundaryLineFactory, _mask);
                }

                // set up
                // ======

                int NoOfEqTotal, NoOfStokesEq, NoOfGaußEq;
                DivergenceFreeBasis DivFreeBasis = this.UseGauß ? new DivergenceFreeBasis(tracker.GridDat, this.Kref, order + 1) : null;
                Basis ScalarBasis = this.UseStokes ? new Basis(this.tracker.GridDat, order + 1) : null; // we loose a order of 1 for the volume rule due to the divergence operator

                NoOfGaußEq   = DivFreeBasis != null ? (DivFreeBasis.Count / D) : 0;
                NoOfStokesEq = (ScalarBasis != null ? ScalarBasis.Length : 0) * D;

                NoOfEqTotal = NoOfGaußEq + NoOfStokesEq;

                // define Nodes
                // ============
                NodeSet NodeSet = null;
                {
                    if (this.Kref.GetType() == typeof(Square))
                    {
                        int K = (int)Math.Ceiling(Math.Sqrt(NoOfEqTotal * 1.7)) + 1;

                        var Nodes1D = GenericBlas.Linspace(-1, 1, K);

                        var _NodeSet = MultidimensionalArray.Create(K * K, 2);
                        int n        = 0;
                        for (int i = 0; i < K /*&& n <= NoOfEq*1.1*/; i++)
                        {
                            for (int j = 0; j < K /*&& n <= NoOfEq*1.1*/; j++)
                            {
                                _NodeSet[n, 0] = Nodes1D[i];
                                _NodeSet[n, 1] = Nodes1D[j];
                                n++;
                            }
                        }

                        NodeSet = new NodeSet(this.Kref, _NodeSet);
                    }
                    else
                    {
                        for (int o = 1; o < 1000000; o++)
                        {
                            var qr = Kref.GetBruteForceQuadRule(o, 0);
                            if (qr.NoOfNodes >= (NoOfEqTotal * 1.1))
                            {
                                NodeSet = qr.Nodes;
                                break;
                            }
                        }
                    }
                }
                int NoOfNodes = NodeSet.GetLength(0);

                // find RHS integrals
                // ==================

                MultidimensionalArray RHS_Gauß = null;
                if (this.UseGauß)
                {
                    stpwGetQuadRuleSet_GaussRHS.Start();
                    RHS_Gauß = this.GaußAnsatzRHS(DivFreeBasis, cellBndSchme, _mask, order);
                    stpwGetQuadRuleSet_GaussRHS.Stop();
                    Debug.Assert(RHS_Gauß.Dimension == 2);
                    Debug.Assert(RHS_Gauß.GetLength(0) == NoOfGaußEq);
                    Debug.Assert(RHS_Gauß.GetLength(1) == _mask.NoOfItemsLocally);
                }

                MultidimensionalArray RHS_Stokes = null;
                if (this.UseStokes)
                {
                    stpwGetQuadRuleSet_StokesRHS.Start();
                    RHS_Stokes = this.StokesAnsatzRHS(ScalarBasis, cellBndLineSchme, _mask, order);
                    stpwGetQuadRuleSet_StokesRHS.Stop();
                    Debug.Assert(RHS_Stokes.Dimension == 2);
                    Debug.Assert(RHS_Stokes.GetLength(0) == NoOfStokesEq);
                    Debug.Assert(RHS_Stokes.GetLength(1) == _mask.NoOfItemsLocally);
                }

                // construct da rule!
                // ==================
                ChunkRulePair <QuadRule>[] SurfaceRule;
                {
                    SurfaceRule = new ChunkRulePair <QuadRule> [_mask.NoOfItemsLocally];
                    var grddat = this.tracker.GridDat;

                    // loop over cells in subgrid...
                    int jSub = 0;
                    foreach (int jCell in _mask.ItemEnum)  // loop over cells in the mask


                    // setup System
                    // ============

                    {
                        NodeSet surfNodes;
                        if (this.SurfaceNodesOnZeroLevset)
                        {
                            surfNodes = ProjectOntoLevset(jCell, NodeSet);
                        }
                        else
                        {
                            surfNodes = NodeSet;
                        }

                        MultidimensionalArray metrics;
                        {
                            int iKref = grddat.Cells.GetRefElementIndex(jCell);
                            metrics = LevelSetData.GetLevelSetNormalReferenceToPhysicalMetrics(surfNodes, jCell, 1);
                        }

                        MultidimensionalArray Mtx_Gauss = null;
                        if (this.UseGauß)
                        {
                            Mtx_Gauss = GaußAnsatzMatrix(DivFreeBasis, surfNodes, jCell);
                            Debug.Assert(Mtx_Gauss.Dimension == 2);
                            Debug.Assert(Mtx_Gauss.GetLength(0) == NoOfGaußEq);
                            Debug.Assert(Mtx_Gauss.GetLength(1) == NoOfNodes);
                        }

                        MultidimensionalArray Mtx_Stokes = null;
                        if (this.UseStokes)
                        {
                            Mtx_Stokes = this.StokesAnsatzMatrix(ScalarBasis, surfNodes, jCell);
                            Debug.Assert(Mtx_Stokes.Dimension == 2);
                            Debug.Assert(Mtx_Stokes.GetLength(0) == NoOfStokesEq);
                            Debug.Assert(Mtx_Stokes.GetLength(1) == NoOfNodes);
                            for (int i = 0; i < NoOfStokesEq; i++)
                            {
                                for (int j = 0; j < NoOfNodes; j++)
                                {
                                    Mtx_Stokes[i, j] /= metrics[0, j];
                                }
                            }
                        }

                        stpwGetQuadRuleSet_SolveRHS.Start();

                        // convert to FORTRAN order
                        MultidimensionalArray _Mtx = MultidimensionalArray.Create(NoOfEqTotal, NoOfNodes);
                        if (this.UseGauß)
                        {
                            _Mtx.ExtractSubArrayShallow(new int[] { 0, 0 }, new int[] { NoOfGaußEq - 1, NoOfNodes - 1 }).Set(Mtx_Gauss);
                        }
                        if (this.UseStokes)
                        {
                            _Mtx.ExtractSubArrayShallow(new int[] { NoOfGaußEq, 0 }, new int[] { NoOfStokesEq + NoOfGaußEq - 1, NoOfNodes - 1 }).Set(Mtx_Stokes);
                        }


                        // convert to FORTRAN order
                        Debug.Assert(NoOfNodes >= NoOfEqTotal);
                        MultidimensionalArray _RHS = MultidimensionalArray.Create(NoOfNodes, 1); // this is also output, so it must be larger!
                        {
                            for (int i = 0; i < NoOfGaußEq; i++)
                            {
                                _RHS[i, 0] = RHS_Gauß[i, jSub];
                            }
                            for (int i = 0; i < NoOfStokesEq; i++)
                            {
                                _RHS[i + NoOfGaußEq, 0] = RHS_Stokes[i, jSub];
                            }
                        }

                        MultidimensionalArray __RHS = null, __Mtx = null;
                        if (this.Docheck)
                        {
                            // values used for testing:
                            __RHS = MultidimensionalArray.Create(NoOfEqTotal);
                            for (int i = 0; i < NoOfEqTotal; i++)
                            {
                                __RHS[i] = _RHS[i, 0];
                            }
                            __Mtx = MultidimensionalArray.Create(_Mtx.NoOfRows, _Mtx.NoOfCols);
                            __Mtx.SetMatrix(_Mtx);
                        }

                        // solve system
                        // ============

                        //int M = _Mtx.NoOfRows;
                        //int N = _Mtx.NoOfCols;
                        //LAPACK.F77_LAPACK.DGELSY(M, N, _Mtx.Entries, _RHS.Entries, 1, 1.0e-14);
                        _Mtx.LeastSquareSolve(_RHS);


                        if (this.Docheck)
                        {
                            // Probe:
                            MultidimensionalArray X = MultidimensionalArray.Create(NoOfNodes); // weights
                            X.ResizeShallow(NoOfNodes, 1).SetMatrix(_RHS);
                            __RHS.Multiply(-1.0, __Mtx, X, 1.0, "j", "jk", "k");
                            double L2_ERR = __RHS.L2Norm();
                            if (L2_ERR > 1.0e-7)
                            {
                                throw new ApplicationException("Quadrature rule in cell " + jCell + " seems to be not very precise: L2_ERR = " + L2_ERR);
                            }
                            //Debug.Assert(L2_ERR < 1.0e-8, "Quadrature rule in cell " + jCell + " seems to be not very precise: L2_ERR = " + L2_ERR);
                            //if (L2_ERR > 1.0e-9)
                            //    Console.WriteLine("Warning: Quadrature rule in cell " + jCell + ": L2_ERR = " + L2_ERR);
                        }


                        stpwGetQuadRuleSet_SolveRHS.Stop();

                        // return da rule!
                        // ===============

                        {
                            {
                                // the surface rule
                                // ----------------

                                QuadRule qr_l = new QuadRule()
                                {
                                    OrderOfPrecision = order,
                                    Weights          = MultidimensionalArray.Create(NoOfNodes),
                                    Nodes            = surfNodes
                                };


                                for (int k = 0; k < NoOfNodes; k++)
                                {
                                    qr_l.Weights[k] = _RHS[k, 0] / metrics[0, k];
                                }

                                SurfaceRule[jSub] = new ChunkRulePair <QuadRule>(Chunk.GetSingleElementChunk(jCell), qr_l);
                            }
                        }
                        jSub++;
                    }
                }


                stpwGetQuadRuleSet.Stop();

                if (!this.m_SurfaceRules.ContainsKey(original_order))
                {
                    this.m_SurfaceRules.Add(original_order, SurfaceRule);
                }
            }
        }
Beispiel #5
0
        /// <summary>
        /// Constructs suitable quadrature rules cells in
        /// <paramref name="mask"/>.
        /// </summary>
        /// <param name="mask">
        /// Cells for which quadrature rules shall be created
        /// </param>
        /// <param name="order">
        /// Desired order of the moment-fitting system. Assuming that
        /// <see cref="surfaceRuleFactory"/> integrates the basis polynomials
        /// exactly over the zero iso-contour (which it usually
        /// doesn't!), the resulting quadrature rules will be exact up to this
        /// order.
        /// </param>
        /// <returns>A set of quadrature rules</returns>
        /// <remarks>
        /// Since the selected level set is generally discontinuous across cell
        /// boundaries, this method does not make use of the fact that
        /// neighboring cells share edges. That is, the optimization will be
        /// performed twice for each inner edge in <paramref name="mask"/>.
        /// </remarks>
        public IEnumerable <IChunkRulePair <QuadRule> > GetQuadRuleSet(ExecutionMask mask, int order)
        {
            using (var tr = new FuncTrace()) {
                CellMask cellMask = mask as CellMask;
                if (cellMask == null)
                {
                    throw new ArgumentException("Mask must be a volume mask", "mask");
                }

                // Note: This is a parallel call, so do this early to avoid parallel confusion
                localCellIndex2SubgridIndex = new SubGrid(cellMask).LocalCellIndex2SubgridIndex;

                int maxLambdaDegree = order + 1;
                int noOfLambdas     = GetNumberOfLambdas(maxLambdaDegree);
                int noOfEdges       = LevelSetData.GridDat.Grid.RefElements[0].NoOfFaces;
                int D = RefElement.SpatialDimension;

                // Get the basis polynomials and integrate them analytically
                Polynomial[] basePolynomials = RefElement.GetOrthonormalPolynomials(order).ToArray();
                Polynomial[] polynomials     = new Polynomial[basePolynomials.Length * D];
                for (int i = 0; i < basePolynomials.Length; i++)
                {
                    Polynomial p = basePolynomials[i];

                    for (int d = 0; d < D; d++)
                    {
                        Polynomial pNew = p.CloneAs();
                        for (int j = 0; j < p.Coeff.Length; j++)
                        {
                            pNew.Exponents[j, d]++;
                            pNew.Coeff[j] /= pNew.Exponents[j, d];
                            pNew.Coeff[j] /= D; // Make sure divergence is Phi again
                        }
                        polynomials[i * D + d] = pNew;
                    }
                }

                // basePolynomials[i] == div(polynomials[i*D], ... , polynomials[i*D + D - 1])
                lambdaBasis = new PolynomialList(polynomials);


                if (RestrictNodes)
                {
                    trafos = new AffineTrafo[mask.NoOfItemsLocally];

                    foreach (Chunk chunk in mask)
                    {
                        foreach (var cell in chunk.Elements.AsSmartEnumerable())
                        {
                            CellMask singleElementMask = new CellMask(
                                LevelSetData.GridDat, Chunk.GetSingleElementChunk(cell.Value));

                            LineAndPointQuadratureFactory.LineQRF lineFactory = this.edgeRuleFactory as LineAndPointQuadratureFactory.LineQRF;
                            if (lineFactory == null)
                            {
                                throw new Exception();
                            }
                            var lineRule  = lineFactory.GetQuadRuleSet(singleElementMask, order).Single().Rule;
                            var pointRule = lineFactory.m_Owner.GetPointFactory().GetQuadRuleSet(singleElementMask, order).Single().Rule;

                            // Also add point rule points since line rule points
                            // are constructed from Gauss rules that do not include
                            // the end points
                            BoundingBox box = new BoundingBox(lineRule.Nodes);
                            box.AddPoints(pointRule.Nodes);

                            int noOfRoots = pointRule.Nodes.GetLength(0);
                            if (noOfRoots <= 1)
                            {
                                // Cell is considered cut because the level set
                                // is very close, but actually isn't. Note that
                                // we can NOT omit the cell (as in the surface
                                // case) as it will be missing in the list of
                                // uncut cells, i.e. this cell would be ignored
                                // completely
                                trafos[localCellIndex2SubgridIndex[cell.Value]] =
                                    AffineTrafo.Identity(RefElement.SpatialDimension);
                                continue;
                            }
                            else if (noOfRoots == 2)
                            {
                                // Go a bit into the direction of the normal
                                // from the center between the nodes in order
                                // not to miss regions with strong curvature
                                double[] center = box.Min.CloneAs();
                                center.AccV(1.0, box.Max);
                                center.ScaleV(0.5);
                                NodeSet centerNode = new NodeSet(RefElement, center);
                                centerNode.LockForever();

                                MultidimensionalArray normal = LevelSetData.GetLevelSetReferenceNormals(centerNode, cell.Value, 1);
                                MultidimensionalArray dist   = LevelSetData.GetLevSetValues(centerNode, cell.Value, 1);

                                double scaling = Math.Sqrt(LevelSetData.GridDat.Cells.JacobiDet[cell.Value]);

                                double[] newPoint = new double[D];
                                for (int d = 0; d < D; d++)
                                {
                                    newPoint[d] = center[d] - normal[0, 0, d] * dist[0, 0] / scaling;
                                }

                                box.AddPoint(newPoint);

                                // Make sure points stay in box
                                for (int d = 0; d < D; d++)
                                {
                                    box.Min[d] = Math.Max(box.Min[d], -1);
                                    box.Max[d] = Math.Min(box.Max[d], 1);
                                }
                            }

                            MultidimensionalArray preImage = RefElement.Vertices.ExtractSubArrayShallow(
                                new int[] { 0, 0 }, new int[] { D, D - 1 });

                            MultidimensionalArray image = MultidimensionalArray.Create(D + 1, D);
                            image[0, 0] = box.Min[0]; // Top left
                            image[0, 1] = box.Max[1];
                            image[1, 0] = box.Max[0]; // Top right
                            image[1, 1] = box.Max[1];
                            image[2, 0] = box.Min[0]; // Bottom left;
                            image[2, 1] = box.Min[1];

                            AffineTrafo trafo = AffineTrafo.FromPoints(preImage, image);
                            trafos[localCellIndex2SubgridIndex[cell.Value]] = trafo;
                        }
                    }
                }

                LambdaCellBoundaryQuadrature cellBoundaryQuadrature =
                    new LambdaCellBoundaryQuadrature(this, edgeRuleFactory, cellMask);
                cellBoundaryQuadrature.Execute();

                LambdaLevelSetSurfaceQuadrature surfaceQuadrature =
                    new LambdaLevelSetSurfaceQuadrature(this, surfaceRuleFactory, cellMask);
                surfaceQuadrature.Execute();

                // Must happen _after_ all parallel calls (e.g., definition of
                // the sub-grid or quadrature) in order to avoid problems in
                // parallel runs
                if (mask.NoOfItemsLocally == 0)
                {
                    var empty = new ChunkRulePair <QuadRule> [0];
                    return(empty);
                }

                if (cachedRules.ContainsKey(order))
                {
                    order = cachedRules.Keys.Where(cachedOrder => cachedOrder >= order).Min();
                    CellMask cachedMask = new CellMask(mask.GridData, cachedRules[order].Select(p => p.Chunk).ToArray());

                    if (cachedMask.Equals(mask))
                    {
                        return(cachedRules[order]);
                    }
                    else
                    {
                        throw new NotImplementedException(
                                  "Case not yet covered yet in combination with caching; deactivate caching to get rid of this message");
                    }
                }

                double[,] quadResults = cellBoundaryQuadrature.Results;
                foreach (Chunk chunk in mask)
                {
                    for (int i = 0; i < chunk.Len; i++)
                    {
                        int iSubGrid = localCellIndex2SubgridIndex[chunk.i0 + i];

                        switch (jumpType)
                        {
                        case JumpTypes.Heaviside:
                            for (int k = 0; k < noOfLambdas; k++)
                            {
                                quadResults[iSubGrid, k] -= surfaceQuadrature.Results[iSubGrid, k];
                            }
                            break;

                        case JumpTypes.OneMinusHeaviside:
                            for (int k = 0; k < noOfLambdas; k++)
                            {
                                quadResults[iSubGrid, k] += surfaceQuadrature.Results[iSubGrid, k];
                            }
                            break;

                        case JumpTypes.Sign:
                            for (int k = 0; k < noOfLambdas; k++)
                            {
                                quadResults[iSubGrid, k] -= 2.0 * surfaceQuadrature.Results[iSubGrid, k];
                            }
                            break;

                        default:
                            throw new NotImplementedException();
                        }
                    }
                }

                BitArray voidCellsArray = new BitArray(LevelSetData.GridDat.Cells.NoOfLocalUpdatedCells);
                BitArray fullCellsArray = new BitArray(LevelSetData.GridDat.Cells.NoOfLocalUpdatedCells);
                foreach (Chunk chunk in cellMask)
                {
                    foreach (var cell in chunk.Elements)
                    {
                        double rhsL2Norm = 0.0;
                        for (int k = 0; k < noOfLambdas; k++)
                        {
                            double entry = quadResults[localCellIndex2SubgridIndex[cell], k];
                            rhsL2Norm += entry * entry;
                        }

                        if (rhsL2Norm < 1e-14)
                        {
                            // All integrals are zero => cell not really cut
                            // (level set is tangent) and fully in void region
                            voidCellsArray[cell] = true;
                            continue;
                        }

                        double l2NormFirstIntegral = quadResults[localCellIndex2SubgridIndex[cell], 0];
                        l2NormFirstIntegral *= l2NormFirstIntegral;
                        double rhsL2NormWithoutFirst = rhsL2Norm - l2NormFirstIntegral;

                        // Beware: This check is only sensible if basis is orthonormal on RefElement!
                        if (rhsL2NormWithoutFirst < 1e-14 &&
                            Math.Abs(l2NormFirstIntegral - RefElement.Volume) < 1e-14)
                        {
                            // All integrals are zero except integral over first integrand
                            // If basis is orthonormal, this implies that cell is uncut and
                            // fully in non-void region since then
                            // \int_K \Phi_i dV = \int_A \Phi_i dV = \delta_{0,i}
                            // However, we have to compare RefElement.Volume since
                            // integration is performed in reference coordinates!
                            fullCellsArray[cell] = true;
                        }
                    }
                }

                var result = new List <ChunkRulePair <QuadRule> >(cellMask.NoOfItemsLocally);

                CellMask emptyCells = new CellMask(LevelSetData.GridDat, voidCellsArray);
                foreach (Chunk chunk in emptyCells)
                {
                    foreach (int cell in chunk.Elements)
                    {
                        QuadRule emptyRule = QuadRule.CreateEmpty(RefElement, 1, RefElement.SpatialDimension);
                        emptyRule.Nodes.LockForever();
                        result.Add(new ChunkRulePair <QuadRule>(
                                       Chunk.GetSingleElementChunk(cell), emptyRule));
                    }
                }

                CellMask fullCells = new CellMask(LevelSetData.GridDat, fullCellsArray);
                foreach (Chunk chunk in fullCells)
                {
                    foreach (int cell in chunk.Elements)
                    {
                        QuadRule fullRule = RefElement.GetQuadratureRule(order);
                        result.Add(new ChunkRulePair <QuadRule>(
                                       Chunk.GetSingleElementChunk(cell), fullRule));
                    }
                }

                CellMask realCutCells = cellMask.Except(emptyCells).Except(fullCells);
                if (RestrictNodes)
                {
                    foreach (Chunk chunk in realCutCells)
                    {
                        foreach (int cell in chunk.Elements)
                        {
                            CellMask singleElementMask = new CellMask(
                                LevelSetData.GridDat, Chunk.GetSingleElementChunk(cell));

                            AffineTrafo trafo = trafos[localCellIndex2SubgridIndex[cell]];
                            Debug.Assert(Math.Abs(trafo.Matrix.Determinant()) > 1e-10);

                            NodeSet nodes       = GetNodes(noOfLambdas).CloneAs();
                            NodeSet mappedNodes = new NodeSet(RefElement, trafo.Transform(nodes));
                            mappedNodes.LockForever();

                            // Remove nodes in negative part
                            MultidimensionalArray levelSetValues  = LevelSetData.GetLevSetValues(mappedNodes, cell, 1);
                            List <int>            nodesToBeCopied = new List <int>(mappedNodes.GetLength(0));
                            for (int n = 0; n < nodes.GetLength(0); n++)
                            {
                                if (levelSetValues[0, n] >= 0.0)
                                {
                                    nodesToBeCopied.Add(n);
                                }
                            }

                            NodeSet reducedNodes = new NodeSet(
                                this.RefElement, nodesToBeCopied.Count, D);
                            for (int n = 0; n < nodesToBeCopied.Count; n++)
                            {
                                for (int d = 0; d < D; d++)
                                {
                                    reducedNodes[n, d] = mappedNodes[nodesToBeCopied[n], d];
                                }
                            }
                            reducedNodes.LockForever();

                            QuadRule optimizedRule = GetOptimizedRule(
                                cell,
                                trafo,
                                reducedNodes,
                                quadResults,
                                order);

                            result.Add(new ChunkRulePair <QuadRule>(
                                           singleElementMask.Single(), optimizedRule));
                        }
                    }
                }
                else
                {
                    // Use same nodes in all cells
                    QuadRule[] optimizedRules = GetOptimizedRules(
                        realCutCells, GetNodes(noOfLambdas), quadResults, order);
                    int ruleIndex = 0;
                    foreach (Chunk chunk in realCutCells)
                    {
                        foreach (var cell in chunk.Elements)
                        {
                            result.Add(new ChunkRulePair <QuadRule>(
                                           Chunk.GetSingleElementChunk(cell), optimizedRules[ruleIndex]));
                            ruleIndex++;
                        }
                    }
                }

                cachedRules[order] = result.OrderBy(p => p.Chunk.i0).ToArray();
                return(cachedRules[order]);
            }
        }