static double JumpNorm(DGField f) { GridData grd = (GridData)f.GridDat; int D = grd.SpatialDimension; var e2cTrafo = grd.Edges.Edge2CellTrafos; f.MPIExchange(); double Unorm = 0; EdgeQuadrature.GetQuadrature( new int[] { D + 1 }, grd, (new EdgeQuadratureScheme()).Compile(grd, f.Basis.Degree * 2), delegate(int i0, int Length, QuadRule QR, MultidimensionalArray EvalResult) { // Evaluate NodeSet NS = QR.Nodes; EvalResult.Clear(); int NoOfNodes = NS.NoOfNodes; for (int j = 0; j < Length; j++) { int iEdge = j + i0; int jCell_IN = grd.Edges.CellIndices[iEdge, 0]; int jCell_OT = grd.Edges.CellIndices[iEdge, 1]; var uDiff = EvalResult.ExtractSubArrayShallow(new int[] { j, 0, 0 }, new int[] { j, NoOfNodes - 1, -1 }); if (jCell_OT >= 0) { int iTrafo_IN = grd.Edges.Edge2CellTrafoIndex[iEdge, 0]; int iTrafo_OT = grd.Edges.Edge2CellTrafoIndex[iEdge, 1]; MultidimensionalArray uIN = MultidimensionalArray.Create(1, NoOfNodes); MultidimensionalArray uOT = MultidimensionalArray.Create(1, NoOfNodes); NodeSet NS_IN = NS.GetVolumeNodeSet(grd, iTrafo_IN); NodeSet NS_OT = NS.GetVolumeNodeSet(grd, iTrafo_OT); f.Evaluate(jCell_IN, 1, NS_IN, uIN); f.Evaluate(jCell_OT, 1, NS_OT, uOT); uDiff.Acc(+1.0, uIN); uDiff.Acc(-1.0, uOT); } else { uDiff.Clear(); } } EvalResult.ApplyAll(x => x * x); }, delegate(int i0, int Length, MultidimensionalArray ResultsOfIntegration) { // SaveIntegrationResults Unorm += ResultsOfIntegration.Sum(); }).Execute(); Unorm = Unorm.MPISum(); return(Unorm.Sqrt()); }
protected override void CorrectLevelSet(int timestepNo, CellMask CellsToCorrect = null, bool ReIteration = false) { DGField Interface_old = Interface.CloneAs(); ScalarFunctionEx FrontTracking_LevelSetDistanceCorrectionFull = delegate(int cell0, int Len, NodeSet Ns, MultidimensionalArray result) { Interface_old.Evaluate(cell0, Len, Ns, result); MultidimensionalArray GlobalNodes = this.Grid.GlobalNodes.GetValue_Cell(Ns, cell0, Len); int c = 0; for (int cell = cell0; cell < cell0 + Len; cell++, c++) { if (NarrowBandCells.Contains(cell)) //if(true) { for (int i = 0; i < Ns.NoOfNodes; i++) { MultidimensionalArray GlobalNode = MultidimensionalArray.Create(1, SpatialDimension); for (int dim = 0; dim < SpatialDimension; dim++) { GlobalNode [0, dim] = GlobalNodes [c, i, dim]; } double Phi_abs = AllEdges[0].SignedDistanceToPoint(GlobalNode); foreach (Edge Edge in AllEdges) { double dist = Edge.SignedDistanceToPoint(GlobalNode); if (Math.Abs(Phi_abs) > Math.Abs(dist)) { Phi_abs = dist; } } result[c, i] = Phi_abs; } } else { for (int i = 0; i < Ns.NoOfNodes; i++) { if (result[c, i] > 0.0) { result[c, i] = 1; } else { result[c, i] = -1; } } } } }; Interface.ProjectField(FrontTracking_LevelSetDistanceCorrectionFull); }
public static ScalarFunctionEx GetCurvatureEnergyFunc(LevelSetTracker LsTrk, DGField Curvature, double sigma, ConventionalDGField[] uI, bool ExtVel, bool squared) { int D = LsTrk.GridDat.SpatialDimension; ScalarFunctionEx CurvatureEnergyFunc = delegate(int i0, int Length, NodeSet nds, MultidimensionalArray result) { Curvature.Evaluate(i0, Length, nds, result); int K = result.GetLength(1); // No nof Nodes MultidimensionalArray U_res = MultidimensionalArray.Create(Length, K, D); for (int i = 0; i < D; i++) { uI.ElementAt(i).Evaluate(i0, Length, nds, U_res.ExtractSubArrayShallow(-1, -1, i)); } var Normals = LsTrk.DataHistories[0].Current.GetLevelSetNormals(nds, i0, Length); for (int j = 0; j < Length; j++) { for (int k = 0; k < K; k++) { double acc = result[j, k]; for (int d = 0; d < D; d++) { // U * N if (!ExtVel) { acc *= U_res[j, k, d] * Normals[j, k, d]; } else { acc *= U_res[j, k, d]; } } if (squared) { result[j, k] = (sigma * acc).Pow2(); } else { result[j, k] = sigma * acc; } } } }; return(CurvatureEnergyFunc); }
protected override void Evaluate(int i0, int Length, QuadRule qr, MultidimensionalArray EvalResult) { if (field is XDGField) { SpeciesId speciesB = new SpeciesId() { cntnt = 11112 }; MultidimensionalArray result = EvalResult.ExtractSubArrayShallow(-1, -1, 0); ((XDGField)field).GetSpeciesShadowField(speciesB).Evaluate(i0, Length, qr.Nodes, result, 0, 0.0); } else { field.Evaluate(i0, Length, qr.Nodes, EvalResult.ExtractSubArrayShallow(-1, -1, 0)); } }
/// <summary> /// Evaluation at a singe point /// </summary> public double Evaluate(Vector Point, DGField F) { var gdat = CellLoc.GrdDat; if (!object.ReferenceEquals(gdat, F.GridDat)) { throw new ArgumentException("Grid mismatch."); } int D = gdat.SpatialDimension; if (Point.Dim != D) { throw new ArgumentException("Spatial dimension mismatch."); } MultidimensionalArray _Point = MultidimensionalArray.Create(1, Point.Dim); _Point.SetRowPt(0, Point); int[] cellIdx = new int[1]; CellLoc.LocalizePointsWithinGrid(_Point, cellIdx, out int NoOfUnassi); if (NoOfUnassi > 0) { throw new ArithmeticException("unable to locate point " + Point); } int j = cellIdx[0]; NodeSet ns = new NodeSet(gdat.iGeomCells.GetRefElement(j), 1, D); bool[] NewtonConvergence = new bool[1]; gdat.TransformGlobal2Local(_Point, ns, j, NewtonConvergence); if (NewtonConvergence[0] == false) { throw new ArithmeticException("unable to transform point " + Point + " to cell " + j + "; Newton method did not converged."); } ns.LockForever(); var R = MultidimensionalArray.Create(1, 1); F.Evaluate(j, 1, ns, R); return(R[0, 0]); }
/// <summary> /// Correction of LevelSet with escaped particles. /// First all particle are assigned with a recalculated level set phi /// Then they are marked as escaped depending on their position in the level set DGField. /// Finally a spherical local level set function is constructed around every escaped particle. /// These functions are overlapped using the following rules: /// For the plus(minus) particles the maximum(minimum) of the function is used. /// Positiv and negative particles are overlapped by prioritizing the level set corrections that are nearer to the interface. /// Method has to be called after the advection step but before the radius adjustment /// </summary> protected override void CorrectLevelSet(int timestepNo, CellMask CellsToCorrect = null, bool ReIteration = false) { DGField Interface_old = Interface.CloneAs(); ScalarFunctionEx Particle_LevelSetCorrection = delegate(int cell0, int Len, NodeSet Ns, MultidimensionalArray result) { Interface_old.Evaluate(cell0, Len, Ns, result); MultidimensionalArray GlobalNodes = this.Grid.GlobalNodes.GetValue_Cell(Ns, cell0, Len); int c = 0; for (int cell = cell0; cell < cell0 + Len; cell++, c++) { if (NarrowBandCells.Contains(cell)) { List <int> listindex = new List <int>(); if (CelltoArrayindex.ContainsKey(cell)) { listindex.Add(cell); } Grid.GetCellNeighbours(cell, GetCellNeighbours_Mode.ViaVertices, out int[] NeighbourIndex, out int[] ConnectingEntities); foreach (int j in NeighbourIndex) { if (CelltoArrayindex.ContainsKey(j)) { listindex.Add(j); } } if (listindex.IsNullOrEmpty()) { continue; } for (int l = 0; l < listindex.Count; l++) { listindex[l] = CelltoArrayindex[listindex[l]]; } for (int i = 0; i < Ns.NoOfNodes; i++) { Dictionary <int, double> Phi_correction = new Dictionary <int, double> { [-1] = result[c, i], [1] = result[c, i] }; if (MinimalDistanceSearch == MinimalDistanceSearchMode.FullSearch) { foreach (int index in listindex) { foreach (KeyValuePair <int, List <SingleLvlSetParticle> > DictEntry in Points[index].ItemList) { foreach (SingleLvlSetParticle SingleParticle in DictEntry.Value) { if (SingleParticle.Escaped == false || SingleParticle.Active == false) { continue; } double dist = 0; for (int dim = 0; dim < SpatialDimension; dim++) { dist += (GlobalNodes[c, i, dim] - SingleParticle.Coordinates[0, dim]).Pow2(); } dist = Math.Sqrt(dist); double Phi_p = DictEntry.Key * (SingleParticle.Radius - dist); Phi_correction[DictEntry.Key] = (DictEntry.Key * Phi_correction[DictEntry.Key] > DictEntry.Key * Phi_p) ? Phi_correction[DictEntry.Key] : Phi_p; } } } } if (Math.Abs(Phi_correction[1]) <= Math.Abs(Phi_correction[-1])) { result[c, i] = Phi_correction[1]; } else { result[c, i] = Phi_correction[-1]; } } } else { for (int i = 0; i < Ns.NoOfNodes; i++) { if (result[c, i] > 0.0) { result[c, i] = 1; } else { result[c, i] = -1; } } } } }; Interface.ProjectField(Particle_LevelSetCorrection); }
public void LimitFieldValues(IEnumerable <DGField> ConservativeVariables, IEnumerable <DGField> DerivedFields) { //IProgram<CNSControl> program; //CNSFieldSet fieldSet = program.WorkingSet; DGField Density = ConservativeVariables.Single(f => f.Identification == Variables.Density.Name); DGField Momentum_0 = ConservativeVariables.Single(f => f.Identification == Variables.Momentum[0].Name); DGField Momentum_1 = ConservativeVariables.Single(f => f.Identification == Variables.Momentum[1].Name); DGField Energy = ConservativeVariables.Single(f => f.Identification == Variables.Energy.Name); foreach (var chunkRulePair in quadRuleSet) { if (chunkRulePair.Chunk.Len > 1) { throw new System.Exception(); } MultidimensionalArray densityValues = MultidimensionalArray.Create(chunkRulePair.Chunk.Len, chunkRulePair.Rule.NoOfNodes); Density.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, densityValues); MultidimensionalArray m0Values = MultidimensionalArray.Create(chunkRulePair.Chunk.Len, chunkRulePair.Rule.NoOfNodes); Momentum_0.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, m0Values); MultidimensionalArray m1Values = MultidimensionalArray.Create(chunkRulePair.Chunk.Len, chunkRulePair.Rule.NoOfNodes); Momentum_1.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, m1Values); MultidimensionalArray energyValues = MultidimensionalArray.Create(chunkRulePair.Chunk.Len, chunkRulePair.Rule.NoOfNodes); Energy.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, energyValues); for (int i = 0; i < chunkRulePair.Chunk.Len; i++) { int cell = i + chunkRulePair.Chunk.i0; CellMask singleCellMask = new CellMask(speciesMap.Tracker.GridDat, chunkRulePair.Chunk); CellQuadratureScheme singleCellScheme = new CellQuadratureScheme( new FixedRuleFactory <QuadRule>(chunkRulePair.Rule), singleCellMask); var singleCellRule = singleCellScheme.Compile(speciesMap.Tracker.GridDat, 99); SpeciesId species = speciesMap.Tracker.GetSpeciesId(speciesMap.Control.FluidSpeciesName); double volume = speciesMap.QuadSchemeHelper.NonAgglomeratedMetrics.CutCellVolumes[species][cell]; double integralDensity = Density.LxError((ScalarFunction)null, (X, a, b) => a, singleCellRule); double meanDensity = integralDensity / volume; if (meanDensity < epsilon) { throw new System.Exception(); } double minDensity = double.MaxValue; for (int j = 0; j < chunkRulePair.Rule.NoOfNodes; j++) { minDensity = Math.Min(minDensity, densityValues[i, j]); } double thetaDensity = Math.Min(1.0, (meanDensity - epsilon) / (meanDensity - minDensity)); if (thetaDensity < 1.0) { Console.WriteLine("Scaled density in cell {0}!", cell); } // Scale for positive density (Beware: Basis is not orthonormal on sub-volume!) for (int j = 0; j < Density.Basis.Length; j++) { Density.Coordinates[cell, j] *= thetaDensity; } Density.AccConstant(meanDensity * (1.0 - thetaDensity), singleCellMask); // Re-evaluate since inner energy has changed densityValues.Clear(); Density.Evaluate(cell, 1, chunkRulePair.Rule.Nodes, densityValues); #if DEBUG // Probe 1 for (int j = 0; j < chunkRulePair.Rule.NoOfNodes; j++) { if (densityValues[i, j] - epsilon < -1e-15) { throw new System.Exception(); } } #endif double integralMomentumX = Momentum_0.LxError((ScalarFunction)null, (X, a, b) => a, singleCellRule); double meanMomentumX = integralMomentumX / volume; double integralMomentumY = Momentum_1.LxError((ScalarFunction)null, (X, a, b) => a, singleCellRule); double meanMomentumY = integralMomentumY / volume; double integralEnergy = Energy.LxError((ScalarFunction)null, (X, a, b) => a, singleCellRule); double meanEnergy = integralEnergy / volume; double meanInnerEnergy = meanEnergy - 0.5 * (meanMomentumX * meanMomentumX + meanMomentumY * meanMomentumY) / meanDensity; if (meanInnerEnergy < epsilon) { throw new System.Exception(); } double minInnerEnergy = double.MaxValue; for (int j = 0; j < chunkRulePair.Rule.NoOfNodes; j++) { double innerEnergy = energyValues[i, j] - 0.5 * (m0Values[i, j] * m0Values[i, j] + m1Values[i, j] * m1Values[i, j]) / densityValues[i, j]; minInnerEnergy = Math.Min(minInnerEnergy, innerEnergy); } // Scale for positive inner energy (Beware: Basis is not orthonormal on sub-volume!) double thetaInnerEnergy = Math.Min(1.0, (meanInnerEnergy - epsilon) / (meanInnerEnergy - minInnerEnergy)); if (thetaInnerEnergy < 1.0) { Console.WriteLine("Scaled inner energy in cell {0}!", cell); } for (int j = 0; j < Density.Basis.Length; j++) { Density.Coordinates[chunkRulePair.Chunk.i0 + i, j] *= thetaInnerEnergy; Momentum_0.Coordinates[chunkRulePair.Chunk.i0 + i, j] *= thetaInnerEnergy; Momentum_1.Coordinates[chunkRulePair.Chunk.i0 + i, j] *= thetaInnerEnergy; Energy.Coordinates[chunkRulePair.Chunk.i0 + i, j] *= thetaInnerEnergy; } Density.AccConstant(meanDensity * (1.0 - thetaInnerEnergy), singleCellMask); Momentum_0.AccConstant(meanMomentumX * (1.0 - thetaInnerEnergy), singleCellMask); Momentum_1.AccConstant(meanMomentumY * (1.0 - thetaInnerEnergy), singleCellMask); Energy.AccConstant(meanEnergy * (1.0 - thetaInnerEnergy), singleCellMask); #if DEBUG // Probe 2 densityValues.Clear(); Density.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, densityValues); m0Values.Clear(); Momentum_0.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, m0Values); m1Values.Clear(); Momentum_1.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, m1Values); energyValues.Clear(); Energy.Evaluate(chunkRulePair.Chunk.i0, chunkRulePair.Chunk.Len, chunkRulePair.Rule.Nodes, energyValues); for (int j = 0; j < chunkRulePair.Rule.NoOfNodes; j++) { StateVector state = new StateVector( speciesMap.GetMaterial(1.0), densityValues[i, j], new Vector(m0Values[i, j], m1Values[i, j], 0.0), energyValues[i, j]); if (!state.IsValid) { throw new System.Exception(); } } #endif } } }
static MultidimensionalArray NodalPatchRecovery(int jCell, double[] Nodes, DGField Phi, out AffineTrafo Trafilein) { GridData GridData = (GridData)Phi.GridDat; int[] Neighbours, Edges; GridData.GetCellNeighbours(jCell, GetCellNeighbours_Mode.ViaVertices, out Neighbours, out Edges); if (Neighbours.Length != 8) { throw new NotSupportedException(); } int[] JS = Neighbours; jCell.AddToArray(ref JS); MultidimensionalArray CellCoordinates = MultidimensionalArray.Create(JS.Length, 4, 2); // 1st index: jcell, top, bottom, left, right var Kref = GridData.iGeomCells.RefElements[0]; if (GridData.iGeomCells.RefElements.Length != 1 || Kref.GetType() != typeof(BoSSS.Foundation.Grid.RefElements.Square)) { throw new NotImplementedException(); } for (int f = 0; f < JS.Length; f++) { GridData.TransformLocal2Global(Kref.Vertices, CellCoordinates.ExtractSubArrayShallow(f, -1, -1), JS[f]); } double xMin = CellCoordinates.ExtractSubArrayShallow(-1, -1, 0).Min(); double xMax = CellCoordinates.ExtractSubArrayShallow(-1, -1, 0).Max(); double yMin = CellCoordinates.ExtractSubArrayShallow(-1, -1, 1).Min(); double yMax = CellCoordinates.ExtractSubArrayShallow(-1, -1, 1).Max(); int N = Nodes.Length; double[] xNodes = new double[N]; double[] yNodes = new double[N]; for (int i = 0; i < N; i++) { xNodes[i] = (xMax - xMin) * 0.5 * (Nodes[i] + 1) + xMin; yNodes[i] = (yMax - yMin) * 0.5 * (Nodes[i] + 1) + yMin; } var TrChey2Glob = new AffineTrafo(2); TrChey2Glob.Matrix[0, 0] = (xMax - xMin) * 0.5; TrChey2Glob.Affine[0] = (xMax - xMin) * 0.5 + xMin; TrChey2Glob.Matrix[1, 1] = (yMax - yMin) * 0.5; TrChey2Glob.Affine[1] = (yMax - yMin) * 0.5 + yMin; MultidimensionalArray[] NodeSets_global = new MultidimensionalArray[JS.Length]; NodeSet[] NodeSets = new NodeSet[JS.Length]; List <double> xNodesCell = new List <double>(); List <double> yNodesCell = new List <double>(); int[][] MapTo_i = new int[JS.Length][]; int[][] MapTo_j = new int[JS.Length][]; for (int f = 0; f < JS.Length; f++) { MultidimensionalArray CellCoords = CellCoordinates.ExtractSubArrayShallow(f, -1, -1); double xMinCell = CellCoords.ExtractSubArrayShallow(-1, 0).Min(); double xMaxCell = CellCoords.ExtractSubArrayShallow(-1, 0).Max(); double yMinCell = CellCoords.ExtractSubArrayShallow(-1, 1).Min(); double yMaxCell = CellCoords.ExtractSubArrayShallow(-1, 1).Max(); int iMin = xNodes.FirstIndexWhere(x => x >= xMinCell); int iMax = xNodes.LastIndexWhere(x => x < xMaxCell); int jMin = yNodes.FirstIndexWhere(y => y >= yMinCell); int jMax = yNodes.LastIndexWhere(y => y < yMaxCell); int NxCell = iMax - iMin + 1; int NyCell = jMax - jMin + 1; int K = NxCell * NyCell; NodeSets_global[f] = MultidimensionalArray.Create(K, 2); MapTo_i[f] = new int[K]; MapTo_j[f] = new int[K]; int k = 0; for (int i = iMin; i <= iMax; i++) { for (int j = jMin; j <= jMax; j++) { MapTo_i[f][k] = i; MapTo_j[f][k] = j; NodeSets_global[f][k, 0] = xNodes[i]; NodeSets_global[f][k, 1] = yNodes[j]; Debug.Assert(GridData.Cells.IsInCell(new double[] { xNodes[i], yNodes[j] }, JS[f], null)); k++; } } NodeSets[f] = new NodeSet(GridData.Cells.GetRefElement(jCell), NxCell * NyCell, 2); GridData.TransformGlobal2Local(NodeSets_global[f], NodeSets[f], JS[f], null); NodeSets[f].LockForever(); } var TrGl2Loc = AffineTrafo.FromPoints(NodeSets_global.Last().ExtractSubArrayShallow(new int[] { 0, 0 }, new int[] { 2, 1 }), NodeSets.Last().ExtractSubArrayShallow(new int[] { 0, 0 }, new int[] { 2, 1 })); Trafilein = (TrGl2Loc * TrChey2Glob).Invert(); var Return = MultidimensionalArray.Create(N, N); for (int f = 0; f < JS.Length; f++) { int K = NodeSets[f].GetLength(0); var Result = MultidimensionalArray.Create(1, K); Phi.Evaluate(JS[f], 1, NodeSets[f], Result); for (int k = 0; k < K; k++) { //double x_i = xNodes[MapTo_i[f][k]]; //double y_j = yNodes[MapTo_j[f][k]]; //double val = (x_i*0.3).Pow2() - (y_j*0.3).Pow2(); Debug.Assert(Return[MapTo_i[f][k], MapTo_j[f][k]] == 0.0); Return[MapTo_i[f][k], MapTo_j[f][k]] = Result[0, k]; } } // ---------------------------------- return(Return); }
protected override void Evaluate(int i0, int Len, QuadRule QR, MultidimensionalArray EvalResult) { NodeSet QuadNodes = QR.Nodes; int D = gridData.SpatialDimension; int NoOfNodes = QuadNodes.NoOfNodes; int GAMMA = m_CodomainMap.NoOfVariables; // GAMMA: number of codom variables // Evaluate Domain & Parameter fields // -------------------------------- Field_Eval.Start(); for (int i = 0; i < m_DomainAndParamFields.Length; i++) { if (m_ValueRequired[i]) { DGField _Field = m_DomainAndParamFields[i]; if (_Field != null) { if (_Field is XDGField) { // jump in parameter i at level-set: separate evaluation for both sides var _xField = _Field as XDGField; _xField.GetSpeciesShadowField(this.SpeciesA).Evaluate(i0, Len, QuadNodes, m_FieldValuesNeg[i]); _xField.GetSpeciesShadowField(this.SpeciesB).Evaluate(i0, Len, QuadNodes, m_FieldValuesPos[i]); } else { // no jump at level set: positive and negative limit of parameter i are equal _Field.Evaluate(i0, Len, QuadNodes, m_FieldValuesPos[i]); m_FieldValuesNeg[i].Set(m_FieldValuesPos[i]); } } else { m_FieldValuesPos[i].Clear(); m_FieldValuesNeg[i].Clear(); } } if (m_GradientRequired[i]) { DGField _Field = m_DomainAndParamFields[i]; if (_Field != null) { if (_Field is XDGField) { // jump in parameter i at level-set: separate evaluation for both sides var _xField = _Field as XDGField; _xField.GetSpeciesShadowField(this.SpeciesA).EvaluateGradient(i0, Len, QuadNodes, m_FieldGradientValuesNeg[i]); _xField.GetSpeciesShadowField(this.SpeciesB).EvaluateGradient(i0, Len, QuadNodes, m_FieldGradientValuesPos[i]); } else { // no jump at level set: positive and negative limit of parameter i are equal _Field.EvaluateGradient(i0, Len, QuadNodes, m_FieldGradientValuesPos[i]); m_FieldGradientValuesNeg[i].Set(m_FieldGradientValuesPos[i]); } } else { m_FieldGradientValuesPos[i].Clear(); m_FieldGradientValuesNeg[i].Clear(); } } } Field_Eval.Stop(); // Evaluate level sets and normals // ------------------------------- ParametersAndNormals.Start(); var NoOfLevSets = m_lsTrk.LevelSets.Count; MultidimensionalArray Normals = m_lsTrk.DataHistories[m_LevSetIdx].Current.GetLevelSetNormals(QuadNodes, i0, Len); MultidimensionalArray NodesGlobal = gridData.GlobalNodes.GetValue_Cell(QuadNodes, i0, Len); ParametersAndNormals.Stop(); // Evaluate basis and test functions // --------------------------------- bool[] ReqV = new bool[GAMMA]; bool[] ReqGradV = new bool[GAMMA]; for (int gamma = 0; gamma < GAMMA; gamma++) { if (Koeff_V[gamma] != null) { ReqV[gamma] = true; } if (Koeff_GradV[gamma] != null) { ReqGradV[gamma] = true; } } MultidimensionalArray[] TestValues; // index: codom variable/test function MultidimensionalArray[] TestGradientValues; // index: codom variable/test function int[,] sectionsTest; Basis_Eval.Start(); LECQuadratureLevelSet <IMutableMatrixEx, double[]> .EvalBasis(i0, Len, this.m_CodomainMap.BasisS, ReqV, ReqGradV, out TestValues, out TestGradientValues, out sectionsTest, QuadNodes); Basis_Eval.Stop(); // Evaluate Integral components // ---------------------------- // loop over codomain variables ... for (int gamma = 0; gamma < GAMMA; gamma++) { // prepare parameters // - - - - - - - - - // set Normal's LevSetIntParams _inParams = new LevSetIntParams(); _inParams.Normal = Normals; // set Nodes Global _inParams.X = NodesGlobal; _inParams.time = this.time; _inParams.LsTrk = this.m_lsTrk; _inParams.i0 = i0; Debug.Assert(_inParams.Len == Len); // clear summation buffers // - - - - - - - - - - - - if (Koeff_V[gamma] != null) { Koeff_V[gamma].Clear(); } if (Koeff_GradV[gamma] != null) { Koeff_GradV[gamma].Clear(); } // Evaluate Bilin. forms // - - - - - - - - - - - { EvalComponent(_inParams, gamma, this.m_NonlinLsForm_V[gamma], this.m_NonlinLsForm_V_Watches[gamma], Koeff_V[gamma], m_FieldValuesPos, m_FieldValuesNeg, m_FieldGradientValuesPos, m_FieldGradientValuesNeg, DELTA, Flux_Eval, delegate(INonlinLevelSetForm_V _comp, LevSetIntParams inp, MultidimensionalArray[] uA, MultidimensionalArray[] uB, MultidimensionalArray[] Grad_uA, MultidimensionalArray[] Grad_uB, MultidimensionalArray SumBuf) { _comp.LevelSetForm_V(_inParams, uA, uB, Grad_uA, Grad_uB, SumBuf); }); } { EvalComponent(_inParams, gamma, this.m_NonlinLsForm_GradV[gamma], this.m_NonlinLsForm_GradV_Watches[gamma], Koeff_GradV[gamma], m_FieldValuesPos, m_FieldValuesNeg, m_FieldGradientValuesPos, m_FieldGradientValuesNeg, DELTA, Flux_Eval, delegate(INonlinLevelSetForm_GradV _comp, LevSetIntParams inp, MultidimensionalArray[] uA, MultidimensionalArray[] uB, MultidimensionalArray[] Grad_uA, MultidimensionalArray[] Grad_uB, MultidimensionalArray SumBuf) { _comp.LevelSetForm_GradV(_inParams, uA, uB, Grad_uA, Grad_uB, SumBuf); }); } } // Summation Loops: multiply with test and trial functions // ------------------------------------------------------- int[] offsetCod = new int[GAMMA]; LECQuadratureLevelSet <IMutableMatrixEx, double[]> . CompOffsets(i0, Len, offsetCod, m_CodomainMap); for (int gamma = 0; gamma < GAMMA; gamma++) { // Evaluate Integrand // - - - - - - - - - var TestVal = TestValues[gamma]; var TestGradVal = TestGradientValues[gamma]; int N; if (TestVal != null) { N = TestVal.GetLength(2); } else if (TestGradVal != null) { N = TestGradVal.GetLength(2); } else { N = 0; } Loops.Start(); // affine offset for (int cr = 0; cr < 2; cr++) // loop over negative/positive species { int[] extr0 = new int[] { 0, 0, sectionsTest[gamma, cr] * N + offsetCod[gamma] }; int[] extrE = new int[] { Len - 1, NoOfNodes - 1, extr0[2] + N - 1 }; var SubRes = EvalResult.ExtractSubArrayShallow(extr0, extrE); if (Koeff_V[gamma] != null) { var Sum_Koeff_V_Cr = Koeff_V[gamma].ExtractSubArrayShallow(-1, -1, cr); SubRes.Multiply(1.0, Sum_Koeff_V_Cr, TestVal, 1.0, "jkn", "jk", "jkn"); } if (Koeff_GradV[gamma] != null) { var Sum_Koeff_NablaV_Cr = Koeff_GradV[gamma].ExtractSubArrayShallow(-1, -1, cr, -1); SubRes.Multiply(1.0, Sum_Koeff_NablaV_Cr, TestGradVal, 1.0, "jkn", "jkd", "jknd"); } } Loops.Stop(); } }
/// <summary> /// Defines the force which is integrated over an immersed boundary, /// called by <see cref="IBMQueries.LiftOrDragForce"/> /// </summary> /// <param name="density"></param> /// <param name="momentum"></param> /// <param name="energy"></param> /// <param name="speciesMap"></param> /// <param name="direction">Direction of the force projection, e.g. 0=x-axis, 1=y-axis</param> /// <param name="cutCellMask">Cells intersected by the interface</param> /// <returns></returns> static ScalarFunctionEx GetSurfaceForce( DGField density, VectorField <DGField> momentum, DGField energy, ImmersedSpeciesMap speciesMap, int direction, CellMask cutCellMask) { return(delegate(int j0, int Len, NodeSet nodes, MultidimensionalArray result) { int noOfNodes = nodes.GetLength(0); int D = nodes.SpatialDimension; double Reynolds = speciesMap.Control.ReynoldsNumber; double Mach = speciesMap.Control.MachNumber; double gamma = speciesMap.Control.EquationOfState.HeatCapacityRatio; double MachScaling = gamma * Mach * Mach; MultidimensionalArray rho = MultidimensionalArray.Create(Len, noOfNodes); density.Evaluate(j0, Len, nodes, rho); MultidimensionalArray[] m = new MultidimensionalArray[CNSEnvironment.NumberOfDimensions]; for (int d = 0; d < CNSEnvironment.NumberOfDimensions; d++) { m[d] = MultidimensionalArray.Create(Len, noOfNodes); momentum[d].Evaluate(j0, Len, nodes, m[d]); } MultidimensionalArray rhoE = MultidimensionalArray.Create(Len, noOfNodes); energy.Evaluate(j0, Len, nodes, rhoE); MultidimensionalArray gradRho = MultidimensionalArray.Create(Len, noOfNodes, D); density.EvaluateGradient(j0, Len, nodes, gradRho); MultidimensionalArray gradM = MultidimensionalArray.Create(Len, noOfNodes, D, D); for (int d = 0; d < D; d++) { momentum[d].EvaluateGradient( j0, Len, nodes, gradM.ExtractSubArrayShallow(-1, -1, d, -1), 0, 0.0); } MultidimensionalArray normals = speciesMap.Tracker.DataHistories[0].Current.GetLevelSetNormals(nodes, j0, Len); Vector3D mVec = new Vector3D(); for (int i = 0; i < Len; i++) { for (int j = 0; j < noOfNodes; j++) { for (int d = 0; d < CNSEnvironment.NumberOfDimensions; d++) { mVec[d] = m[d][i, j]; } Material material = speciesMap.GetMaterial(double.NaN); StateVector state = new StateVector(material, rho[i, j], mVec, rhoE[i, j]); double mu = 0.0; if (Reynolds != 0.0) { mu = state.GetViscosity(j0 + j) / Reynolds; } double[,] gradU = new double[D, D]; for (int d1 = 0; d1 < D; d1++) { for (int d2 = 0; d2 < D; d2++) { // Apply chain rule gradU[d1, d2] = (gradM[i, j, d1, d2] - state.Momentum[d1] / state.Density * gradRho[i, j, d2]) / state.Density; } } double divU = gradU[0, 0] + gradU[1, 1]; switch (direction) { // Attention: Changed sign, because normal vector is pointing inwards, not outwards! case 0: // x-Direction result[i, j, 0] = -state.Pressure / MachScaling * normals[i, j, 0] + mu * (2.0 * gradU[0, 0] - 2.0 / 3.0 * divU) * normals[i, j, 0] //tau_11 * n_1 + mu * (gradU[0, 1] + gradU[1, 0]) * normals[i, j, 1]; //tau_12 * n_2 break; case 1: // y-Direction result[i, j, 0] = -state.Pressure / MachScaling * normals[i, j, 1] + mu * (gradU[0, 1] + gradU[1, 0]) * normals[i, j, 0] //tau_12 * n_1 + mu * (2.0 * gradU[1, 1] - 2.0 / 3.0 * divU) * normals[i, j, 1]; //tau_22*n_2 break; default: throw new ArgumentException("Lift and Drag currently only in 2D implemented"); } } } }); }
/// <summary> /// L2 error of some quantity derived from the state vector (e.g., /// entropy) with respect to given reference solution. The quadrature /// is determined from the settings in <see cref="IBMControl"/> /// </summary> /// <param name="quantityOfInterest"></param> /// <param name="referenceSolution"></param> /// <returns></returns> public static Query L2Error(Func <StateVector, double> quantityOfInterest, Func <double[], double, double> referenceSolution) { return(delegate(IApplication <AppControl> app, double time) { IProgram <CNSControl> program = app as IProgram <CNSControl>; if (program == null) { throw new Exception(); } ImmersedSpeciesMap speciesMap = program.SpeciesMap as ImmersedSpeciesMap; IBMControl control = program.Control as IBMControl; if (speciesMap == null || control == null) { throw new Exception( "Query is only valid for immersed boundary runs"); } SpeciesId species = speciesMap.Tracker.GetSpeciesId(control.FluidSpeciesName); int order = control.LevelSetQuadratureOrder; CellQuadratureScheme scheme = speciesMap.QuadSchemeHelper.GetVolumeQuadScheme( species, true, speciesMap.SubGrid.VolumeMask); var composititeRule = scheme.Compile(program.GridData, order); IChunkRulePair <QuadRule>[] chunkRulePairs = composititeRule.ToArray(); DGField density = program.WorkingSet.Density; VectorField <DGField> momentum = program.WorkingSet.Momentum; DGField energy = program.WorkingSet.Energy; // Construct dummy field since L2Error is currently only supported // for Field's; However, _avoid_ a projection. DGField dummy = new SinglePhaseField(new Basis(program.GridData, 0)); Material material = speciesMap.GetMaterial(double.NaN); int index = 0; double value = dummy.LxError( (ScalarFunctionEx) delegate(int j0, int Len, NodeSet nodes, MultidimensionalArray result) { MultidimensionalArray input = program.GridData.GlobalNodes.GetValue_Cell(nodes, j0, Len); Chunk chunk = chunkRulePairs[index].Chunk; QuadRule rule = chunkRulePairs[index].Rule; if (chunk.i0 != j0 || chunk.Len != Len) { throw new Exception(); } if (rule.NoOfNodes != nodes.GetLength(0)) { throw new Exception(); } MultidimensionalArray rho = MultidimensionalArray.Create(chunk.Len, rule.NoOfNodes); density.Evaluate(chunk.i0, chunk.Len, nodes, rho); MultidimensionalArray[] m = new MultidimensionalArray[CNSEnvironment.NumberOfDimensions]; for (int d = 0; d < CNSEnvironment.NumberOfDimensions; d++) { m[d] = MultidimensionalArray.Create(chunk.Len, rule.NoOfNodes); momentum[d].Evaluate(chunk.i0, chunk.Len, nodes, m[d]); } MultidimensionalArray rhoE = MultidimensionalArray.Create(chunk.Len, rule.NoOfNodes); energy.Evaluate(chunk.i0, chunk.Len, nodes, rhoE); double[] X = new double[CNSEnvironment.NumberOfDimensions]; Vector3D mVec = new Vector3D(); for (int i = 0; i < chunk.Len; i++) { for (int j = 0; j < rule.NoOfNodes; j++) { for (int d = 0; d < CNSEnvironment.NumberOfDimensions; d++) { X[d] = input[i, j, d]; mVec[d] = m[d][i, j]; } StateVector state = new StateVector(material, rho[i, j], mVec, rhoE[i, j]); double qoi = quantityOfInterest(state); Debug.Assert( !double.IsNaN(qoi), "Encountered node with unphysical state" + " (not able to determine quantity of interest)"); result[i, j] = qoi - referenceSolution(X, time); } } index++; }, (X, a, b) => (a - b) * (a - b), composititeRule); // No value is NaN, but the results. How can this be? // => All values around 0, but values in void region are a little // farther away from the exact solution // => weights in the void zone sum up to something slightly negative Debug.Assert( value >= 0, "Encountered unphysical norm even though individual values where valid." + " This indicates a problem with cut-cell quadrature."); return Math.Sqrt(value); }); }