private double PerformSurfaceQuadrature(Modes mode, IQuadRuleFactory <QuadRule> volumeFactory, IQuadRuleFactory <QuadRule> edgeFactory, SubGrid cutCellGrid, int order, Stopwatch timer) { using (new FuncTrace()) { CellQuadratureScheme volInstr = new CellQuadratureScheme( volumeFactory, cutCellGrid.VolumeMask); CellBoundaryQuadratureScheme edgeInstr = new CellBoundaryQuadratureScheme( new CellBoundaryFromEdgeRuleFactory <CellBoundaryQuadRule>( GridData, Grid.RefElements[0], edgeFactory), cutCellGrid.VolumeMask); ScalarFieldLevelSetIntegrator quadrature = new ScalarFieldLevelSetIntegrator( levelSetTracker, SinglePhaseField, volInstr.Compile(GridData, order), edgeInstr.Compile(GridData, order), order, cutCellGrid, 0); timer.Start(); double result = quadrature.ExecuteA().Storage.Sum(); timer.Stop(); return(result); } }
MultidimensionalArray StokesAnsatzRHS(Basis TestBasis, CellBoundaryQuadratureScheme cellBndSchme, CellMask _mask, int order) { var GridDat = this.tracker.GridDat; CellBoundaryQuadrature <CellBoundaryQuadRule> qBnd = null; int N = TestBasis.Length; int D = GridDat.SpatialDimension; MultidimensionalArray RHS = MultidimensionalArray.Create(D, N, _mask.NoOfItemsLocally); double[] CellN = new double[D]; // cell normal double[] SurfN = new double[D]; // level-set normal double[] OutwardTang = new double[D]; // level-set tangent, outward of cell if (D != 2) { throw new NotSupportedException("Currently only supported for spatial dimension of 2."); } //MultidimensionalArray Nudes = null; int jSgrd = 0; qBnd = CellBoundaryQuadrature <CellBoundaryQuadRule> .GetQuadrature(new int[] { D, N }, GridDat, cellBndSchme.Compile(GridDat, order), delegate(int i0, int Length, CellBoundaryQuadRule NS, MultidimensionalArray EvalResult) { // Evaluate //MultidimensionalArray BasisValues = TestBasis.Evaluate(0); // reference //var LSNormals = LsTrk.GetLevelSetReferenceNormals(iLevSet, 0, i0, Length); // reference MultidimensionalArray BasisValues = TestBasis.CellEval(NS.Nodes, i0, Length); // physical MultidimensionalArray LSNormals = this.LevelSetData.GetLevelSetNormals(NS.Nodes, i0, Length); // physical for (int i = 0; i < Length; i++) // loop over cells //if(i0 + i == 1) { // EvalResult.ExtractSubArrayShallow(i, -1, -1, -1).Clear(); // continue; //} { CellBoundaryQuadRule cR = qBnd.CurrentRule; int[] NodesPerEdge = cR.NumbersOfNodesPerFace; var Kref = cR.RefElement; int NoOfFaces = Kref.NoOfFaces; int iNode = 0; Debug.Assert(NoOfFaces == NodesPerEdge.Length); for (int e = 0; e < NoOfFaces; e++) // loop over the faces of the cell { if (NodesPerEdge[e] <= 0) { continue; } // reference: //for (int d = 0; d < D; d++) { // CellN[d] = Kref.FaceNormals[e, d]; //} // ~~~~ // physical: var FaceNodes = new NodeSet(Kref, cR.Nodes.ExtractSubArrayShallow(new int[] { iNode, 0 }, new int[] { iNode + NodesPerEdge[e] - 1, D - 1 })); var FaceNormals = MultidimensionalArray.Create(NodesPerEdge[e], D); GridDat.Edges.GetNormalsForCell(FaceNodes, i0, e, FaceNormals); // ~~~~ for (int _n = 0; _n < NodesPerEdge[e]; _n++) // loop over nodes in one edge { for (int d = 0; d < D; d++) { SurfN[d] = LSNormals[i, iNode, d]; CellN[d] = FaceNormals[_n, d]; // physical } tangente(SurfN, CellN, OutwardTang); for (int n = 0; n < N; n++) // loop over Test polynomials (the same as the basis polynomials) { for (int d = 0; d < D; d++) // loop over spatial direction { EvalResult[i, iNode, d, n] = BasisValues[i, iNode, n] * OutwardTang[d]; // physical //EvalResult[i, iNode, d, n] = BasisValues[iNode, n]*OutwardTang[d]; // reference } } iNode++; } } Debug.Assert(iNode == EvalResult.GetLength(1)); } }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { // SaveIntegrationResults for (int i = 0; i < Length; i++) { var ResPart = RHS.ExtractSubArrayShallow(new int[] { 0, 0, jSgrd }, new int[] { D - 1, N - 1, jSgrd - 1 }); int NoOfFaces = ResultsOfIntegration.GetLength(1); for (int e = 0; e < NoOfFaces; e++) { var ip = ResultsOfIntegration.ExtractSubArrayShallow(new int[] { i, e, 0, 0 }, new int[] { i - 1, e - 1, D - 1, N - 1 }); ResPart.Acc(1.0, ip); } jSgrd++; } }, cs : CoordinateSystem.Physical); qBnd.Execute(); var ret = RHS.ResizeShallow(N * D, _mask.NoOfItemsLocally); return(ret); }
protected MultidimensionalArray GaußAnsatzRHS(DivergenceFreeBasis TestBasis, CellBoundaryQuadratureScheme cellBndScheme, CellMask _mask, int order) { var _Context = this.tracker.GridDat; int D = this.tracker.GridDat.Grid.SpatialDimension; int N = TestBasis.Count; var coordSys = CoordinateSystem.Reference; var LsTrk = this.tracker; int Nrhs = _mask.NoOfItemsLocally; Debug.Assert(N % D == 0); N /= D; MultidimensionalArray RHS = MultidimensionalArray.Create(N, Nrhs); var splx = this.Kref; int NoOfFaces = splx.NoOfFaces; //var normals = _Context.GridDat.Normals; CellBoundaryQuadrature <CellBoundaryQuadRule> qBnd = null; int jSgrd = 0; qBnd = CellBoundaryQuadrature <CellBoundaryQuadRule> .GetQuadrature(new int[] { N }, _Context, cellBndScheme.Compile(_Context, order), delegate(int i0, int Length, CellBoundaryQuadRule QR, MultidimensionalArray EvalResult) { // Evaluate NodeSet Nodes = QR.Nodes; MultidimensionalArray BasisValues; if (coordSys == CoordinateSystem.Physical) { //BasisValues = TestBasis.CellEval(Nodes, i0, Length); throw new NotImplementedException("todo"); } else if (coordSys == CoordinateSystem.Reference) { BasisValues = TestBasis.Values.GetValues(Nodes); } else { throw new NotImplementedException(); } for (int i = 0; i < Length; i++) // loop over cells { CellBoundaryQuadRule cR = qBnd.CurrentRule; int[] NodesPerEdge = cR.NumbersOfNodesPerFace; Debug.Assert(object.ReferenceEquals(splx, cR.RefElement)); int iNode = 0; Debug.Assert(NoOfFaces == NodesPerEdge.Length); for (int e = 0; e < NoOfFaces; e++) // loop over the faces of the cell { for (int _n = 0; _n < NodesPerEdge[e]; _n++) // loop over nodes in one edge { for (int n = 0; n < N; n++) // loop over Test polynomials (the same as the basis polynomials) { double acc = 0; for (int d = 0; d < D; d++) // loop over spatial directions { if (coordSys == CoordinateSystem.Physical) { throw new NotImplementedException("todo"); //int q = _Context.GridDat.LocalCellIndexToEdges[i+i0, e]; //int iEdge = Math.Abs(q) - 1; //double Nsign = Math.Sign(q); //double Nd = normals[iEdge, d]; //EvalResult[i, iNode, n, d] = BasisValues[i, iNode, n]*Nd*Nsign; } else { Debug.Assert(coordSys == CoordinateSystem.Reference); double Nd = splx.FaceNormals[e, d]; //Debug.Assert(Nd == normals[iEdge, d]*Nsign); acc += BasisValues[iNode, n *D + d] * Nd; } } EvalResult[i, iNode, n] = acc; } iNode++; } } Debug.Assert(iNode == EvalResult.GetLength(1)); } }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { // SaveIntegrationResults for (int i = 0; i < Length; i++) { var ResPart = RHS.ExtractSubArrayShallow(new int[] { 0, jSgrd }, new int[] { N - 1, jSgrd - 1 }); for (int e = 0; e < NoOfFaces; e++) { var ip = ResultsOfIntegration.ExtractSubArrayShallow(new int[] { i, e, 0 }, new int[] { i - 1, e - 1, N - 1 }); ResPart.Acc(1.0, ip); } jSgrd++; } }, cs : coordSys); qBnd.Execute(); return(RHS); }
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); } } }
/// <summary> /// Compares the cell surface and the boundary integral. /// </summary> public static void TestSealing(IGridData gdat) { int J = gdat.iLogicalCells.NoOfLocalUpdatedCells; int[,] E2Clog = gdat.iGeomEdges.LogicalCellIndices; // // compute cell surface via edge integrals // MultidimensionalArray CellSurf1 = MultidimensionalArray.Create(J); EdgeQuadrature.GetQuadrature(new int[] { 1 }, gdat, new EdgeQuadratureScheme().Compile(gdat, 2), delegate(int i0, int Length, QuadRule rule, MultidimensionalArray EvalResult) { // evaluate EvalResult.SetAll(1.0); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { // save results for (int i = 0; i < Length; i++) { int iEdge = i + i0; int jCellIn = E2Clog[iEdge, 0]; int jCellOt = E2Clog[iEdge, 1]; CellSurf1[jCellIn] += ResultsOfIntegration[i, 0]; if (jCellOt >= 0 && jCellOt < J) { CellSurf1[jCellOt] += ResultsOfIntegration[i, 0]; } } }).Execute(); // // compute cell surface via cell boundary integrals // var cbqs = new CellBoundaryQuadratureScheme(false, null); foreach (var Kref in gdat.iGeomCells.RefElements) { cbqs.AddFactory(new StandardCellBoundaryQuadRuleFactory(Kref)); } MultidimensionalArray CellSurf2 = MultidimensionalArray.Create(J); CellBoundaryQuadrature <CellBoundaryQuadRule> .GetQuadrature(new int[] { 1 }, gdat, cbqs.Compile(gdat, 2), delegate(int i0, int Length, CellBoundaryQuadRule rule, MultidimensionalArray EvalResult) { // evaluate EvalResult.SetAll(1.0); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { // save results int NoOfFaces = ResultsOfIntegration.GetLength(1); for (int i = 0; i < Length; i++) { int jCell = i + i0; for (int iFace = 0; iFace < NoOfFaces; iFace++) { CellSurf2[jCell] += ResultsOfIntegration[i, iFace, 0]; } } }, cs : CoordinateSystem.Physical).Execute(); // // compare // MultidimensionalArray Err = CellSurf1.CloneAs(); Err.Acc(-1.0, CellSurf2); double ErrNorm = Err.L2Norm(); double TotSurf = CellSurf1.L2Norm(); Console.WriteLine("Area Check " + Err.L2Norm()); Assert.LessOrEqual(ErrNorm / TotSurf, 1.0e-8); //for (int j = 0; j < J; j++) { // if (Err[j].Abs() >= 1.0e-6) { // Console.WriteLine("Mismatch between edge area and cell surface area in cell {0}, GlobalId {1}, No. of neighbors {3}, {2:0.####E-00}", j, gdat.CurrentGlobalIdPermutation.Values[j], Err[j], gdat.Cells.CellNeighbours[j].Length); // schas.SetMeanValue(j, 1); // } //} }