void CalculateDerivative() { DuDx.Clear(); DuDy.Clear(); DvDx.Clear(); DvDy.Clear(); DuDx.Derivative(1.0, m_app.WorkingSet.Velocity.Current[0], 0); DuDy.Derivative(1.0, m_app.WorkingSet.Velocity.Current[0], 1); DvDx.Derivative(1.0, m_app.WorkingSet.Velocity.Current[1], 0); DvDy.Derivative(1.0, m_app.WorkingSet.Velocity.Current[1], 1); if (m_app.GridData.SpatialDimension == 3) { DuDz.Clear(); DvDz.Clear(); DwDx.Clear(); DwDy.Clear(); DwDz.Clear(); DuDz.Derivative(1.0, m_app.WorkingSet.Velocity.Current[0], 2); DvDz.Derivative(1.0, m_app.WorkingSet.Velocity.Current[1], 2); DwDx.Derivative(1.0, m_app.WorkingSet.Velocity.Current[2], 0); DwDy.Derivative(1.0, m_app.WorkingSet.Velocity.Current[2], 1); DwDz.Derivative(1.0, m_app.WorkingSet.Velocity.Current[2], 2); } }
/// <summary> /// Calculate A Level-Set field from the Explicit description /// </summary> /// <param name="LevelSet">Target Field</param> /// <param name="LsTrk"></param> public void ProjectToDGLevelSet(SinglePhaseField LevelSet, LevelSetTracker LsTrk = null) { CellMask VolMask; if (LsTrk == null) { VolMask = CellMask.GetFullMask(LevelSet.Basis.GridDat); } else { // set values in positive and negative FAR region to +1 and -1 CellMask Near = LsTrk.Regions.GetNearMask4LevSet(0, 1); CellMask PosFar = LsTrk.Regions.GetLevelSetWing(0, +1).VolumeMask.Except(Near); CellMask NegFar = LsTrk.Regions.GetLevelSetWing(0, -1).VolumeMask.Except(Near); LevelSet.Clear(PosFar); LevelSet.AccConstant(1, PosFar); LevelSet.Clear(NegFar); LevelSet.AccConstant(-1, NegFar); // project Fourier levelSet to DGfield on near field VolMask = Near; } LevelSet.Clear(VolMask); // scalar function is already vectorized for parallel execution // nodes in global coordinates LevelSet.ProjectField(1.0, PhiEvaluation(mode), new Foundation.Quadrature.CellQuadratureScheme(true, VolMask)); // check the projection error projErr_phiDG = LevelSet.L2Error( PhiEvaluation(mode), new Foundation.Quadrature.CellQuadratureScheme(true, VolMask)); //if (projErr_phiDG >= 1e-5) // Console.WriteLine("WARNING: LevelSet projection error onto PhiDG = {0}", projErr_phiDG); // project on higher degree field and take the difference //SinglePhaseField higherLevSet = new SinglePhaseField(new Basis(LevelSet.GridDat, LevelSet.Basis.Degree * 2), "higherLevSet"); //higherLevSet.ProjectField(1.0, PhiEvaluator(), new Foundation.Quadrature.CellQuadratureScheme(true, VolMask)); //double higherProjErr = higherLevSet.L2Error( // PhiEvaluator(), // new Foundation.Quadrature.CellQuadratureScheme(true, VolMask)); //Console.WriteLine("LevelSet projection error onto higherPhiDG = {0}", higherProjErr.ToString()); // check the projection from current sample points on the DGfield //MultidimensionalArray interP = MultidimensionalArray.Create(current_interfaceP.Lengths); //if (this is PolarFourierLevSet) { // interP = ((PolarFourierLevSet)this).interfaceP_cartesian; //} else { // interP = current_interfaceP; //} //projErr_interface = InterfaceProjectionError(LevelSet, current_interfaceP); //if (projErr_interface >= 1e-3) // Console.WriteLine("WARNING: Interface projection Error onto PhiDG = {0}", projErr_interface); }
/// <summary> /// Sets the positive Far-field of the level-set to +1 and the negative side to -1 /// </summary> /// <param name="LevelSet"></param> /// <param name="Domain"></param> /// <param name="PosMask"></param> public void SetFarField(SinglePhaseField LevelSet, CellMask Domain, CellMask PosMask) { // set cells outside narrow band to +/- 1 var Pos = PosMask.Except(Domain); var Neg = PosMask.Complement().Except(Domain); LevelSet.Clear(Pos); LevelSet.AccConstant(1, Pos); LevelSet.Clear(Neg); LevelSet.AccConstant(-1, Neg); }
/// <summary> /// Initialize logging of energy for Orr-Sommerfeld problem /// </summary> void InitLogEnergyOrrSommerfeld() { Energy = new SinglePhaseField(new Basis(base.GridData, 20)); if (base.MPIRank == 0) { Log_Energy = base.DatabaseDriver.FsDriver.GetNewLog("PerturbationEnergy", CurrentSessionInfo.ID); } Energy.Clear(); Energy.ProjectFunction( 1.0, (X, U, cell) => (1.0 - X[1] * X[1] - U[0]) * (1.0 - X[1] * X[1] - U[0]) + U[1] * U[1], null, WorkingSet.Velocity.Current[0], WorkingSet.Velocity.Current[1]); Energy_t0 = Energy.IntegralOver(null); if (base.MPIRank == 0) { Log_Energy.Write("Time\tPerturbationEnergy\tOmega"); Log_Energy.WriteLine(); Log_Energy.Write("{0}\t{1}\t{2}", (0.0).ToString("0.000", NumberFormatInfo.InvariantInfo), Energy_t0.ToString("0.0000000000E+00", NumberFormatInfo.InvariantInfo), (2.234976E-03).ToString("0.000000E+00", NumberFormatInfo.InvariantInfo)); Log_Energy.Flush(); } if (base.MPIRank == 0) { Console.WriteLine("E(t0) = " + Energy_t0); } }
public void MakeContinuous(SinglePhaseField DGLevelSet, SinglePhaseField LevelSet, CellMask Domain) { FEMLevSet.ProjectDGField(1.0, DGLevelSet, Domain); LevelSet.Clear(); FEMLevSet.AccToDGField(1.0, LevelSet, Domain); LevelSet.AccLaidBack(1.0, DGLevelSet, Domain.Complement()); }
/// <summary> /// Rhs corrector, cf. right-hand side of Eq. (17) in /// B. Klein, F. Kummer, M. Keil, and M. Oberlack, /// An extension of the SIMPLE based discontinuous Galerkin solver to unsteady incompressible flows, J. Comput. Phys., 2013. /// </summary> /// <param name="dt"></param> /// <param name="SpatialComponent"></param> /// <returns></returns> protected override IList <double> DefineRhs(double dt, int SpatialComponent) { double[] Rhs = new double[m_DivB4.DOFLocal]; // calculate mass defect m_DivB4.Clear(); SolverUtils.CalculateMassDefect_Divergence(m_VelocityDivergence, m_Velocity_Intrmed, m_DivB4); Rhs.AccV(1.0, m_DivB4.CoordinateVector); // pressure stabilization if (base.m_solverConf.Control.PressureStabilizationScaling > 0.0) { double[] PressureStabi = new double[Rhs.Length]; m_PressureStabilization.OperatorMatrix.SpMVpara(1.0, m_Pressure.CoordinateVector, 0.0, PressureStabi); Rhs.AccV(1.0, PressureStabi); } // reference point pressure if (base.m_solverConf.Control.PressureReferencePoint != null) { BoSSS.Solution.NSECommon.SolverUtils.SetRefPtPressure_Rhs(Rhs, base.m_solverConf.PressureReferencePointIndex, m_MatAsmblyCorrector.i0); } return(Rhs); }
/// <summary> /// Log energy for Orr-Sommerfeld problem. /// </summary> void LogEnergyOrrSommerfeld(int TimestepNo, double phystime, double dt) { Energy.Clear(); Energy.ProjectFunction( 1.0, (X, U, cell) => (1.0 - X[1] * X[1] - U[0]) * (1.0 - X[1] * X[1] - U[0]) + U[1] * U[1], null, WorkingSet.Velocity.Current[0], WorkingSet.Velocity.Current[1]); EnergyIntegral = Energy.IntegralOver(null); omega = 1 / (2 * (phystime + dt)) * Math.Log(EnergyIntegral / Energy_t0); if (base.MPIRank == 0) { Log_Energy.WriteLine(); Log_Energy.Write("{0}\t{1}\t{2}", (phystime + dt).ToString("0.000", NumberFormatInfo.InvariantInfo), EnergyIntegral.ToString("0.0000000000E+00", NumberFormatInfo.InvariantInfo), omega.ToString("0.000000E+00", NumberFormatInfo.InvariantInfo)); Log_Energy.Flush(); if (TimestepNo == base.Control.NoOfTimesteps) { Log_Energy.Close(); } } }
protected override IList <double> DefineRhs(double dt, int SpatialComponent) { double[] Rhs = new double[DivB4.DOFLocal]; // calculate mass defect DivB4.Clear(); SolverUtils.CalculateMassDefect_Divergence(VelocityDivergence, Velocity_Intrmed, DivB4); Rhs.AccV(1.0, DivB4.CoordinateVector); // time derivative density if (base.m_solverConf.Control.Algorithm == SolutionAlgorithms.Unsteady_SIMPLE) { SinglePhaseField DensitySummand = new SinglePhaseField(DivB4.Basis); BDF.ComputeDensitySummand(dt, base.m_solverConf.BDFOrder, Temperature, EoS, DensitySummand); Rhs.AccV(1.0, DensitySummand.CoordinateVector); } // reference point pressure if (base.m_solverConf.Control.PressureReferencePoint != null) { BoSSS.Solution.NSECommon.SolverUtils.SetRefPtPressure_Rhs(Rhs, base.m_solverConf.PressureReferencePointIndex, MatAsmblyCorrector.i0); } return(Rhs); }
/// <summary> /// Sets level-set and solution at time (<paramref name="time"/> + <paramref name="dt"/>). /// </summary> double DelUpdateLevelset(DGField[] CurrentState, double time, double dt, double UnderRelax, bool _incremental) { // new time double t = time + dt; // project new level-set LevSet.ProjectField(this.Control.LevelSet.Vectorize(t)); LsTrk.UpdateTracker(incremental: _incremental); // exact solution for new timestep uEx.GetSpeciesShadowField("A").ProjectField(NonVectorizedScalarFunction.Vectorize(this.Control.uEx_A, t)); uEx.GetSpeciesShadowField("B").ProjectField(NonVectorizedScalarFunction.Vectorize(this.Control.uEx_B, t)); u.Clear(); u.Acc(1.0, uEx); // markieren, wo ueberhaupt A und B sind Amarker.Clear(); Bmarker.Clear(); Amarker.AccConstant(+1.0, LsTrk._Regions.GetSpeciesSubGrid("A").VolumeMask); Bmarker.AccConstant(1.0, LsTrk._Regions.GetSpeciesSubGrid("B").VolumeMask); // MPI rank MPICellRank.Clear(); MPICellRank.AccConstant(base.MPIRank); // return level-set residual return(0.0); }
/// <summary> /// Sets level-set and solution at time (<paramref name="time"/> + <paramref name="dt"/>). /// </summary> double DelUpdateLevelset(DGField[] CurrentState, double time, double dt, double UnderRelax, bool _incremental) { // new time double t = time + dt; // project new level-set double s = 1.0; LevSet.ProjectField((x, y) => - (x - s * t).Pow2() - y.Pow2() + (2.4).Pow2()); LsTrk.UpdateTracker(incremental: _incremental); LsTrk.PushStacks(); // exact solution for new timestep uXEx.GetSpeciesShadowField("A").ProjectField((x, y) => x + alpha_A * t); uXEx.GetSpeciesShadowField("B").ProjectField((x, y) => x + alpha_B * t); // uX.Clear(); //uX.Acc(1.0, uXEx); // single-phase-properties u.Clear(); u.ProjectField(NonVectorizedScalarFunction.Vectorize(uEx, t)); Grad_u.Clear(); Grad_u.Gradient(1.0, u); MagGrad_u.Clear(); MagGrad_u.ProjectFunction(1.0, (double[] X, double[] U, int jCell) => Math.Sqrt(U[0].Pow2() + U[1].Pow2()), new Foundation.Quadrature.CellQuadratureScheme(), Grad_u.ToArray()); // return level-set residual return(0.0); }
private double SetUpConfiguration() { testCase.UpdateLevelSet(levelSet); levelSetTracker.UpdateTracker(__NearRegionWith: 0, incremental: false, __LevSetAllowedMovement: 2); XDGField.Clear(); XDGField.GetSpeciesShadowField("A").ProjectField( 1.0, testCase.JumpingFieldSpeciesAInitialValue, default(CellQuadratureScheme)); XDGField.GetSpeciesShadowField("B").ProjectField( 1.0, testCase.JumpingFieldSpeciesBInitialValue, default(CellQuadratureScheme)); SinglePhaseField.Clear(); SinglePhaseField.ProjectField(testCase.ContinuousFieldInitialValue); double referenceValue = testCase.Solution; if (testCase is IVolumeTestCase) { QuadRule standardRule = Grid.RefElements[0].GetQuadratureRule(2 * XDGField.Basis.Degree + 1); ScalarFieldQuadrature uncutQuadrature = new ScalarFieldQuadrature( GridData, XDGField, new CellQuadratureScheme( new FixedRuleFactory <QuadRule>(standardRule), levelSetTracker.Regions.GetCutCellSubGrid().Complement().VolumeMask), standardRule.OrderOfPrecision); uncutQuadrature.Execute(); referenceValue -= uncutQuadrature.Result; } return(referenceValue); }
/// <summary> /// Sets level-set and solution at time (<paramref name="time"/> + <paramref name="dt"/>). /// </summary> double DelUpdateLevelset(DGField[] CurrentState, double time, double dt, double UnderRelax, bool _incremental) { // new time double t = time + dt; // project new level-set double s = 1.0; LevSet.ProjectField((x, y) => - (x - s * t).Pow2() - y.Pow2() + (2.4).Pow2()); LsTrk.UpdateTracker(incremental: _incremental); // exact solution for new timestep uEx.GetSpeciesShadowField("A").ProjectField((x, y) => x + alpha_A * t); uEx.GetSpeciesShadowField("B").ProjectField((x, y) => x + alpha_B * t); u.Clear(); u.Acc(1.0, uEx); // markieren, wo ueberhaupt A und B sind Amarker.Clear(); Bmarker.Clear(); Amarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("A").VolumeMask); Bmarker.AccConstant(1.0, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); // MPI rank MPICellRank.Clear(); MPICellRank.AccConstant(base.MPIRank); // return level-set residual return(0.0); }
/// <summary> /// Determines the total curvature which is twice the mean curvature /// Please pay attention to the sign of this expression! /// </summary> /// <param name="Output">The total curvature</param> /// <param name="optionalSubGrid"> /// Subgrid which can be defined, for example for carrying out /// computations on a narrow band /// </param> /// <param name="bndMode"> /// Definition of the behavior at subgrid boundaries /// </param> public void ComputeTotalCurvatureByFlux(SinglePhaseField Output, SubGrid optionalSubGrid = null, SpatialOperator.SubGridBoundaryModes bndMode = SpatialOperator.SubGridBoundaryModes.OpenBoundary) { if (this.m_Basis.Degree <= 1) { throw new ArgumentException("For correct computation of these level set quantities, the level set has to be at least of degree 2!"); } Basis basisForNormalVec = new Basis(this.GridDat, this.m_Basis.Degree - 1); Basis basisForCurvature = new Basis(this.GridDat, this.m_Basis.Degree - 2); Func <Basis, string, SinglePhaseField> fac = (Basis b, string id) => new SinglePhaseField(b, id); VectorField <SinglePhaseField> normalVector = new VectorField <SinglePhaseField>(this.GridDat.SpatialDimension, basisForNormalVec, fac); ComputeNormalByFlux(normalVector, optionalSubGrid, bndMode); VectorField <SinglePhaseField> secondDerivatives = new VectorField <SinglePhaseField>(this.GridDat.SpatialDimension, basisForCurvature, fac); Output.Clear(); for (int i = 0; i < normalVector.Dim; i++) { secondDerivatives[i].DerivativeByFlux(1.0, normalVector[i], i, optionalSubGrid, bndMode); Output.Acc(-1.0, secondDerivatives[i], optionalSubGrid.VolumeMask); } }
/// <summary> /// Calculate Extension /// </summary> public void ConstructExtension(IList <DGField> InterfaceValue = null, bool nearfield = false) { using (new FuncTrace()) { Extension.Clear(nearfield ? LevelSetTracker.Regions.GetNearFieldMask(1) : null); ComputeMatrices(InterfaceValue ?? InterfaceParams, nearfield); //Solve System double[] RHS = OpAffine.CloneAs(); RHS.ScaleV(-1.0); ISparseSolver slv = Control.solverFactory(); if (nearfield) { SubGrid subGrid = LevelSetTracker.Regions.GetNearFieldSubgrid(1); int[] SubVecIdx = Extension.Mapping.GetSubvectorIndices(subGrid, true, new int[] { 0 }); int L = SubVecIdx.Length; MsrMatrix SubMatrix = new MsrMatrix(L); double[] SubRHS = new double[L]; double[] SubSolution = new double[L]; OpMatrix.AccSubMatrixTo(1.0, SubMatrix, SubVecIdx, default(int[]), SubVecIdx, default(int[])); SubRHS.Clear(); SubSolution.Clear(); SubRHS.AccV(1.0, RHS, default(int[]), SubVecIdx); SubSolution.AccV(1.0, Extension.CoordinateVector, default(int[]), SubVecIdx); slv.DefineMatrix(SubMatrix); slv.Solve(SubSolution, SubRHS); Extension.Clear(subGrid.VolumeMask); Extension.CoordinateVector.AccV(1.0, SubSolution, SubVecIdx, default(int[])); } else { slv.DefineMatrix(OpMatrix); slv.Solve(Extension.CoordinateVector, RHS); } slv.Dispose(); } }
/// <summary> /// see <see cref="DGField.Laplacian(double,DGField,CellMask)"/>; /// </summary> public override void Laplacian(double alpha, DGField f, CellMask em) { SinglePhaseField tmp = new SinglePhaseField(this.Basis, "tmp"); int D = this.GridDat.SpatialDimension; for (int d = 0; d < D; d++) { tmp.Clear(); tmp.Derivative(1.0, f, d, em); this.Derivative(alpha, tmp, d, em); } }
public void MakeContinuous(SinglePhaseField DGLevelSet, SinglePhaseField LevelSet, CellMask Domain) { if (ReferenceEquals(DGLevelSet, LevelSet)) { // Do nothing } else { LevelSet.Clear(); LevelSet.AccLaidBack(1.0, DGLevelSet); } }
//static int PlotCont = 1; void UpdateMarkerFields() { CutMarker.Clear(); NearMarker.Clear(); DOFMarker.Clear(); foreach (int j in this.LsTrk.Regions.GetCutCellMask4LevSet(0).ItemEnum) { CutMarker.SetMeanValue(j, 1); } foreach (int j in this.LsTrk.Regions.GetNearFieldMask(1).ItemEnum) { NearMarker.SetMeanValue(j, 1); } int J = this.GridData.iLogicalCells.NoOfLocalUpdatedCells; for (int j = 0; j < J; j++) { DOFMarker.SetMeanValue(j, this.u.Basis.GetLength(j)); } /* * Tecplot.PlotFields(new DGField[] { CutMarker, NearMarker }, "Markers-" + PlotCont + ".csv", 0.0, 1); * LsTrk.Regions.GetCutCellMask().SaveToTextFile("Cut-" + PlotCont + ".csv", false); * LsTrk.Regions.GetSpeciesMask("A").SaveToTextFile("SpcA-" + PlotCont + ".csv", false); * LsTrk.Regions.GetSpeciesMask("B").SaveToTextFile("SpcB-" + PlotCont + ".csv", false); * * int qOrd = this.LinearQuadratureDegree; * var sch = LsTrk.GetXDGSpaceMetrics(LsTrk.SpeciesIdS.ToArray(), qOrd, 1).XQuadSchemeHelper; * * var schCut = sch.GetLevelSetquadScheme(0, LsTrk.Regions.GetCutCellMask()); * var RuleCut = schCut.SaveCompile(this.GridData, qOrd); * ICompositeQuadRule_Ext.SumOfWeightsToTextFileVolume(RuleCut, this.GridData, "CutRule-" + PlotCont + ".csv"); * * var schB = sch.GetVolumeQuadScheme(LsTrk.GetSpeciesId("B")); * var RuleB = schB.SaveCompile(this.GridData, qOrd); * ICompositeQuadRule_Ext.SumOfWeightsToTextFileVolume(RuleB, this.GridData, "B_Rule-" + PlotCont + ".csv"); * * var schA = sch.GetVolumeQuadScheme(LsTrk.GetSpeciesId("A")); * var RuleA = schA.SaveCompile(this.GridData, qOrd); * ICompositeQuadRule_Ext.SumOfWeightsToTextFileVolume(RuleA, this.GridData, "A_Rule-" + PlotCont + ".csv"); * * var eschB = sch.GetEdgeQuadScheme(LsTrk.GetSpeciesId("B")); * var ERuleB = eschB.SaveCompile(this.GridData, qOrd); * ICompositeQuadRule_Ext.SumOfWeightsToTextFileEdge(ERuleB, this.GridData, "Be_Rule-" + PlotCont + ".csv"); * * var eschA = sch.GetEdgeQuadScheme(LsTrk.GetSpeciesId("A")); * var ERuleA = eschA.SaveCompile(this.GridData, qOrd); * ICompositeQuadRule_Ext.SumOfWeightsToTextFileEdge(ERuleA, this.GridData, "Ae_Rule-" + PlotCont + ".csv"); * * PlotCont++; */ }
void DerivativeByFluxLinear(SinglePhaseField fin, SinglePhaseField fres, int d, SinglePhaseField fBnd) { var Op = (new LinearDerivFlx(d)).Operator(); BlockMsrMatrix OpMtx = new BlockMsrMatrix(fres.Mapping, fin.Mapping); double[] OpAff = new double[fres.Mapping.LocalLength]; Op.ComputeMatrixEx(fin.Mapping, new DGField[] { fBnd }, fres.Mapping, OpMtx, OpAff, OnlyAffine: false); fres.Clear(); fres.CoordinateVector.Acc(1.0, OpAff); OpMtx.SpMV(1.0, fin.CoordinateVector, 1.0, fres.CoordinateVector); }
private void Filter(SinglePhaseField FiltrdField, int NoOfSweeps, CellMask CC) { Basis patchRecoveryBasis = FiltrdField.Basis; L2PatchRecovery l2pr = new L2PatchRecovery(patchRecoveryBasis, patchRecoveryBasis, CC, true); SinglePhaseField F_org = FiltrdField.CloneAs(); for (int pass = 0; pass < NoOfSweeps; pass++) { F_org.Clear(); F_org.Acc(1.0, FiltrdField); FiltrdField.Clear(); l2pr.Perform(FiltrdField, F_org); } }
public static void ProjectArtificalViscosityToDGField(SinglePhaseField avField, PerssonSensor sensor, double sensorLimit, double maxViscosity, CellMask cellMask = null) { MultidimensionalArray h_min = avField.GridDat.iGeomCells.h_min; int p = sensor.fieldToTestRestricted.Basis.Degree + 1; if (cellMask == null) { cellMask = CellMask.GetFullMask(avField.GridDat); } avField.Clear(); foreach (int cell in cellMask.ItemEnum) { double localViscosity = GetViscosity(cell, h_min[cell], p, sensor.GetValue(cell), sensorLimit, maxViscosity); avField.SetMeanValue(cell, localViscosity); } }
/// <summary> /// Calculate current Nusselt number. /// </summary> public void CalculateNusseltNumber() { dTdx.Clear(); //dTdx.Derivative(1.0, Temperature, 0, GridDat.BoundaryCells.VolumeMask); dTdx.DerivativeByFlux(1.0, Temperature, 0); dTdy.Clear(); //dTdy.Derivative(1.0, Temperature, 1, GridDat.BoundaryCells.VolumeMask); dTdy.DerivativeByFlux(1.0, Temperature, 1); for (int bc = 0; bc < Nusselt.Length; bc++) { double LocalNusselt = NusseltIntegrals[bc].Evaluate(); double GlobalNusselt = 0.0; unsafe { MPI.Wrappers.csMPI.Raw.Allreduce((IntPtr)(&LocalNusselt), (IntPtr)(&GlobalNusselt), 1, MPI.Wrappers.csMPI.Raw._DATATYPE.DOUBLE, MPI.Wrappers.csMPI.Raw._OP.SUM, MPI.Wrappers.csMPI.Raw._COMM.WORLD); } Nusselt[bc] = GlobalNusselt; } }
/// <summary> /// calculate the curvature corresponding to the Level-Set /// </summary> /// <param name="Curvature"></param> /// <param name="LevSetGradient"></param> /// <param name="VolMask"></param> public void ProjectToDGCurvature(SinglePhaseField Curvature, out VectorField <SinglePhaseField> LevSetGradient, CellMask VolMask = null) { if (VolMask == null) { VolMask = CellMask.GetFullMask(Curvature.Basis.GridDat); } Curvature.Clear(); Curvature.ProjectField(1.0, CurvEvaluation(mode), new Foundation.Quadrature.CellQuadratureScheme(true, VolMask)); // check the projection error projErr_curv = Curvature.L2Error( CurvEvaluation(mode), new Foundation.Quadrature.CellQuadratureScheme(true, VolMask)); //if (projErr_curv >= 1e-5) // Console.WriteLine("WARNING: Curvature projection error onto PhiDG = {0}", projErr_curv); LevSetGradient = null; }
protected override double RunSolverOneStep(int TimestepNo, double phystime, double dt) { Console.WriteLine(" Timestep # " + TimestepNo + ", phystime = " + phystime); //phystime = 1.8; LsUpdate(phystime); // operator-matrix assemblieren MsrMatrix OperatorMatrix = new MsrMatrix(u.Mapping, u.Mapping); double[] Affine = new double[OperatorMatrix.RowPartitioning.LocalLength]; MultiphaseCellAgglomerator Agg; MassMatrixFactory Mfact; // Agglomerator setup int quadOrder = Op.QuadOrderFunction(new int[] { u.Basis.Degree }, new int[0], new int[] { u.Basis.Degree }); //Agg = new MultiphaseCellAgglomerator(new CutCellMetrics(MomentFittingVariant, quadOrder, LsTrk, ), this.THRESHOLD, false); Agg = LsTrk.GetAgglomerator(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, quadOrder, this.THRESHOLD); Console.WriteLine("Inter-Process agglomeration? " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.InterProcessAgglomeration); if (this.THRESHOLD > 0.01) { TestAgglomeration_Extraploation(Agg); TestAgglomeration_Projection(quadOrder, Agg); } // operator matrix assembly Op.ComputeMatrixEx(LsTrk, u.Mapping, null, u.Mapping, OperatorMatrix, Affine, false, 0.0, true, Agg.CellLengthScales, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(OperatorMatrix, Affine, u.Mapping, u.Mapping); // mass matrix factory Mfact = LsTrk.GetXDGSpaceMetrics(new SpeciesId[] { LsTrk.GetSpeciesId("B") }, quadOrder, 1).MassMatrixFactory;// new MassMatrixFactory(u.Basis, Agg); // Mass matrix/Inverse Mass matrix //var MassInv = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, true, LsTrk.GetSpeciesId("B")); var Mass = Mfact.GetMassMatrix(u.Mapping, new double[] { 1.0 }, false, LsTrk.GetSpeciesId("B")); Agg.ManipulateMatrixAndRHS(Mass, default(double[]), u.Mapping, u.Mapping); var MassInv = Mass.InvertBlocks(OnlyDiagonal: true, Subblocks: true, ignoreEmptyBlocks: true, SymmetricalInversion: false); // test that operator depends only on B-species values double DepTest = LsTrk.Regions.GetSpeciesSubGrid("B").TestMatrixDependency(OperatorMatrix, u.Mapping, u.Mapping); Console.WriteLine("Matrix dependency test: " + DepTest); Assert.LessOrEqual(DepTest, 0.0); // diagnostic output Console.WriteLine("Number of Agglomerations (all species): " + Agg.TotalNumberOfAgglomerations); Console.WriteLine("Number of Agglomerations (species 'B'): " + Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).AggInfo.SourceCells.NoOfItemsLocally.MPISum()); // operator auswerten: double[] x = new double[Affine.Length]; BLAS.daxpy(x.Length, 1.0, Affine, 1, x, 1); OperatorMatrix.SpMVpara(1.0, u.CoordinateVector, 1.0, x); MassInv.SpMV(1.0, x, 0.0, du_dx.CoordinateVector); Agg.GetAgglomerator(LsTrk.GetSpeciesId("B")).Extrapolate(du_dx.Mapping); // markieren, wo ueberhaupt A und B sind Bmarker.AccConstant(1.0, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Amarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("A").VolumeMask); Xmarker.AccConstant(+1.0, LsTrk.Regions.GetSpeciesSubGrid("X").VolumeMask); // compute error ERR.Clear(); ERR.Acc(1.0, du_dx_Exact, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); ERR.Acc(-1.0, du_dx, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double L2Err = ERR.L2Norm(LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); Console.WriteLine("L2 Error: " + L2Err); XERR.Clear(); XERR.GetSpeciesShadowField("B").Acc(1.0, ERR, LsTrk.Regions.GetSpeciesSubGrid("B").VolumeMask); double xL2Err = XERR.L2Norm(); Console.WriteLine("L2 Error (in XDG space): " + xL2Err); // check error if (this.THRESHOLD > 0.01) { // without agglomeration, the error in very tiny cut-cells may be large over the whole cell // However, the error in the XDG-space should be small under all circumstances Assert.LessOrEqual(L2Err, 1.0e-6); } Assert.LessOrEqual(xL2Err, 1.0e-6); bool IsPassed = ((L2Err <= 1.0e-6 || this.THRESHOLD <= 0.01) && xL2Err <= 1.0e-7); if (IsPassed) { Console.WriteLine("Test PASSED"); } else { Console.WriteLine("Test FAILED: check errors."); } // return/Ende base.NoOfTimesteps = 17; //base.NoOfTimesteps = 2; dt = 0.3; return(dt); }
/// <summary> /// Obtaining the time integrated spatial discretization of the reinitialization equation in a narrow band around the zero level set, based on a Godunov's numerical Hamiltonian calculation /// </summary> /// <param name="LS"> The level set function </param> /// <param name="Restriction"> The narrow band around the zero level set </param> /// <param name="NumberOfTimesteps"> /// maximum number of pseudo-timesteps /// </param> /// <param name="thickness"> /// The smoothing width of the signum function. /// This is the main stabilization parameter for re-initialization. /// It should be set to approximately 3 cells. /// </param> /// <param name="TimestepSize"> /// size of the pseudo-timestep /// </param> public void ReInitialize(LevelSet LS, SubGrid Restriction, double thickness, double TimestepSize, int NumberOfTimesteps) { using (var tr = new FuncTrace()) { // log parameters: tr.Info("thickness: " + thickness.ToString(NumberFormatInfo.InvariantInfo)); tr.Info("TimestepSize: " + TimestepSize.ToString(NumberFormatInfo.InvariantInfo)); tr.Info("NumberOfTimesteps: " + NumberOfTimesteps); ExplicitEuler TimeIntegrator; SpatialOperator SO; Func <int[], int[], int[], int> QuadratureOrder = QuadOrderFunc.NonLinear(3); if (m_ctx.SpatialDimension == 2) { SO = new SpatialOperator(1, 5, 1, QuadratureOrder, new string[] { "LS", "LSCGV", "LSDG[0]", "LSUG[0]", "LSDG[1]", "LSUG[1]", "Result" }); SO.EquationComponents["Result"].Add(new GodunovHamiltonian(m_ctx, thickness)); SO.Commit(); TimeIntegrator = new RungeKutta(m_Scheme, SO, new CoordinateMapping(LS), new CoordinateMapping(LSCGV, LSDG[0], LSUG[0], LSDG[1], LSUG[1]), sgrd: Restriction); } else { SO = new SpatialOperator(1, 7, 1, QuadratureOrder, new string[] { "LS", "LSCGV", "LSDG[0]", "LSUG[0]", "LSDG[1]", "LSUG[1]", "LSDG[2]", "LSUG[2]", "Result" }); SO.EquationComponents["Result"].Add(new GodunovHamiltonian(m_ctx, thickness)); SO.Commit(); TimeIntegrator = new RungeKutta(m_Scheme, SO, new CoordinateMapping(LS), new CoordinateMapping(LSCGV, LSDG[0], LSUG[0], LSDG[1], LSUG[1], LSDG[2], LSUG[2]), sgrd: Restriction); } // Calculating the gradients in each sub-stage of a Runge-Kutta integration procedure ExplicitEuler.ChangeRateCallback EvalGradients = delegate(double t1, double t2) { LSUG.Clear(); CalculateLevelSetGradient(LS, LSUG, "Upwind", Restriction); LSDG.Clear(); CalculateLevelSetGradient(LS, LSDG, "Downwind", Restriction); LSCG.Clear(); CalculateLevelSetGradient(LS, LSCG, "Central", Restriction); LSCGV.Clear(); var VolMask = (Restriction != null) ? Restriction.VolumeMask : null; LSCGV.ProjectAbs(1.0, VolMask, LSCG.ToArray()); }; TimeIntegrator.OnBeforeComputeChangeRate += EvalGradients; { EvalGradients(0, 0); var GodunovResi = new SinglePhaseField(LS.Basis, "Residual"); SO.Evaluate(1.0, 0.0, LS.Mapping, TimeIntegrator.ParameterMapping.Fields, GodunovResi.Mapping, Restriction); //Tecplot.Tecplot.PlotFields(ArrayTools.Cat<DGField>( LSUG, LSDG, LS, GodunovResi), "Residual", 0, 3); } // pseudo-timestepping // =================== double factor = 1.0; double time = 0; LevelSet prevLevSet = new LevelSet(LS.Basis, "prevLevSet"); CellMask RestrictionMask = (Restriction == null) ? null : Restriction.VolumeMask; for (int i = 0; (i < NumberOfTimesteps); i++) { tr.Info("Level set reinitialization pseudo-timestepping, timestep " + i); // backup old Levelset // ------------------- prevLevSet.Clear(); prevLevSet.Acc(1.0, LS, RestrictionMask); // time integration // ---------------- double dt = TimestepSize * factor; tr.Info("dt = " + dt.ToString(NumberFormatInfo.InvariantInfo) + " (factor = " + factor.ToString(NumberFormatInfo.InvariantInfo) + ")"); TimeIntegrator.Perform(dt); time += dt; // change norm // ------ prevLevSet.Acc(-1.0, LS, RestrictionMask); double ChangeNorm = prevLevSet.L2Norm(RestrictionMask); Console.WriteLine("Reinit: PseudoTime: {0} - Changenorm: {1}", i, ChangeNorm); //Tecplot.Tecplot.PlotFields(new SinglePhaseField[] { LS }, m_ctx, "Reinit-" + i, "Reinit-" + i, i, 3); } //*/ } }
/// <summary> /// Single run of the solver /// </summary> protected override double RunSolverOneStep(int TimestepNo, double phystime, double dt) { using (new FuncTrace()) { if (Control.ExactSolution_provided) { Tex.Clear(); Tex.ProjectField(this.Control.InitialValues_Evaluators["Tex"]); RHS.Clear(); RHS.ProjectField(this.Control.InitialValues_Evaluators["RHS"]); } if (Control.AdaptiveMeshRefinement == false) { base.NoOfTimesteps = -1; if (TimestepNo > 1) { throw new ApplicationException("steady-state-equation."); } base.TerminationKey = true; } // Update matrices // --------------- UpdateMatrices(); // call solver // ----------- double mintime, maxtime; bool converged; int NoOfIterations; ClassicSolve(out mintime, out maxtime, out converged, out NoOfIterations); Console.WriteLine("finished; " + NoOfIterations + " iterations."); Console.WriteLine("converged? " + converged); Console.WriteLine("Timespan: " + mintime + " to " + maxtime + " seconds"); base.QueryHandler.ValueQuery("minSolRunT", mintime, true); base.QueryHandler.ValueQuery("maxSolRunT", maxtime, true); base.QueryHandler.ValueQuery("Conv", converged ? 1.0 : 0.0, true); base.QueryHandler.ValueQuery("NoIter", NoOfIterations, true); base.QueryHandler.ValueQuery("NoOfCells", this.GridData.CellPartitioning.TotalLength, true); base.QueryHandler.ValueQuery("DOFs", T.Mapping.TotalLength, true); base.QueryHandler.ValueQuery("BlockSize", T.Basis.Length, true); if (base.Control.ExactSolution_provided) { Error.Clear(); Error.AccLaidBack(1.0, Tex); Error.AccLaidBack(-1.0, T); double L2_ERR = Error.L2Norm(); Console.WriteLine("\t\tL2 error on " + this.Grid.NumberOfCells + ": " + L2_ERR); base.QueryHandler.ValueQuery("SolL2err", L2_ERR, true); } // evaluate residual // ================= { this.ResiualKP1.Clear(); if (this.Control.InitialValues_Evaluators.ContainsKey("RHS")) { this.ResiualKP1.ProjectField(this.Control.InitialValues_Evaluators["RHS"]); } var ev = this.LapaceIp.GetEvaluator(T.Mapping, ResiualKP1.Mapping); ev.Evaluate(-1.0, 1.0, ResiualKP1.CoordinateVector); } // return // ====== return(0.0); } }
/// <summary> /// Based on the Ideas by /// C. Basting and D. Kuzmin, /// “A minimization-based finite element formulation for interface-preserving level set reinitialization”, /// Computing, vol. 95, no. 1, pp. 13–25, Dec. 2012. /// Create Spatial Operators and build the corresponding Matrices /// For the Left-Hand Side of the ReInitProblem /// RHS is computed on the fly in <see cref="ReInitSingleStep"/> /// The Bulk component is constant unless the grid changes, thus it is computed in <see cref="BuildOperators(CellQuadratureScheme)"/>. /// The Interface component changes with its motion. /// This component is calculated in <see cref="UpdateOperators(CellQuadratureScheme)"/>. /// </summary> /// <param name="LSTrck"></param> /// <param name="Control">various parameters <paramref name="EllipticReinitControl"/></param> /// <param name="HMFOrder">order of tghe interface quadrature</param> public EllipticReInit(LevelSetTracker LSTrck, EllipticReInitAlgoControl Control, SinglePhaseField LevelSetForReInit = null) { this.Control = Control; this.LevelSetTracker = LSTrck; if (LevelSetForReInit == null) { Phi = LevelSetTracker.LevelSets[0] as SinglePhaseField; } else { Phi = LevelSetForReInit; } this.underrelaxation = Control.underrelaxation; Residual = new SinglePhaseField(Phi.Basis); OldPhi = new SinglePhaseField(Phi.Basis); NewPhi = new SinglePhaseField(Phi.Basis); foreach (SinglePhaseField f in new List <SinglePhaseField> { Residual, OldPhi, NewPhi }) { f.Clear(); f.Acc(1.0, Phi); } this.D = LevelSetTracker.GridDat.SpatialDimension; this.ConvergenceCriterion = Control.ConvergenceCriterion; this.MaxIteration = Control.MaxIt; double PenaltyBase = ((double)((Phi.Basis.Degree + 1) * (Phi.Basis.Degree + D))) / ((double)D); // Choose Forms according to Upwinding or Central Fluxes string[] paramNames; int noOfParamFields; IEquationComponent BulkForm; RHSForm myRHSForm; LevelSetGradient = new VectorField <SinglePhaseField>(D, Phi.Basis, "LevelSetGradient", SinglePhaseField.Factory); MeanLevelSetGradient = new VectorField <SinglePhaseField>(D, new Basis(Phi.GridDat, 0), "MeanLevelSetGradient", SinglePhaseField.Factory); if (Control.Upwinding) { paramNames = new string[] { "OldLevelSet", "MeanLevelSetGradient[0]", "MeanLevelSetGradient[1]" }; noOfParamFields = D; LevelSetGradient.Clear(); LevelSetGradient.Gradient(1.0, Phi); //LevelSetGradient.GradientByFlux(1.0, Phi); MeanLevelSetGradient.Clear(); MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient); parameterFields = ArrayTools.Cat(new SinglePhaseField[] { OldPhi }, MeanLevelSetGradient.ToArray()); //throw new NotImplementedException("ToDO"); BulkForm = new EllipticReInitUpwindForm_Laplace(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck); myRHSForm = new EllipticReInitUpwindForm_RHS(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck); OldDirection = new double[MeanLevelSetGradient.CoordinateVector.ToArray().Length]; for (int i = 0; i < MeanLevelSetGradient.CoordinateVector.Length; i++) { OldDirection[i] = Math.Sign(MeanLevelSetGradient.CoordinateVector[i]); } NewDirection = OldDirection.CloneAs(); } else { paramNames = new string[] { }; noOfParamFields = 0; parameterFields = new SinglePhaseField[] { }; BulkForm = new CentralDifferencesLHSForm(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck.GridDat.Cells.cj); myRHSForm = new CentralDifferencesRHSForm(Control.PenaltyMultiplierFlux * PenaltyBase, LSTrck); } // SIP for the bulk Phase //this.Operator_bulk = new SpatialOperator(1, noOfParamFields, 1, QuadOrderFunc.SumOfMaxDegrees(1, RoundUp: false), variableNames); this.Operator_bulk = BulkForm.Operator(); // Zero at the Interface // Calculate Quadrature Order Func <int[], int[], int[], int> InterfaceQuadOrder; InterfaceQuadOrder = QuadOrderFunc.FixedOrder(Phi.Basis.Degree * 2 + 2); // Generate Interface Operator this.Operator_interface = (new EllipticReInitInterfaceForm(Control.PenaltyMultiplierInterface * PenaltyBase, LSTrck)).XOperator(new[] { "A" }, InterfaceQuadOrder); // Nonlinear Part on the RHS // switch for the potential functions switch (Control.Potential) { case ReInitPotential.BastingDoubleWell: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.DoubleWell(d, b)); break; }; case ReInitPotential.BastingSingleWell: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.SingleWell(d, b)); break; }; case ReInitPotential.SingleWellNear: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.SingleWellNear(d, b)); break; }; case ReInitPotential.P4DoubleWell: { Console.WriteLine("Warning - This Option for Elliptic ReInit does not work well"); myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.DoubleWellAlternative(d, b)); break; }; case ReInitPotential.SingleWellOnCutDoubleWellElse: { myRHSForm.DiffusionRate = ((d, b) => DiffusionRates.SingleWellOnCutDoubleWellElse(d, b)); break; } } Operator_RHS = myRHSForm.Operator(QuadOrderFunc.SumOfMaxDegrees(2, RoundUp: true)); // The result of the nonlinear part on the rhs is projected on a single-phase field RHSField = new SinglePhaseField(Phi.Basis, "RHS"); OpMatrix = new MsrMatrix(this.Phi.Mapping, this.Phi.Mapping); OpAffine = new double[OpMatrix.RowPartitioning.LocalLength]; // Matrix and RHS for the Bulk component OpMatrix_bulk = new MsrMatrix(this.Phi.Mapping, this.Phi.Mapping); OpAffine_bulk = new double[OpMatrix.RowPartitioning.LocalLength]; // Matrix and RHS for the Interface Penalty OpMatrix_interface = new MsrMatrix(this.Phi.Mapping, this.Phi.Mapping); OpAffine_interface = new double[OpMatrix.RowPartitioning.LocalLength]; // Init Parameter Fields OldPhi.Clear(); OldPhi.Acc(1.0, Phi); // Compute Matrices UpdateBulkMatrix(); }
/// <summary> /// One Iteration of the ReInitialization /// Operators must be built first /// </summary> /// <param name="ChangeRate"> /// L2-Norm of the Change-Rate in the level set in this reinit step /// </param> /// <param name="Restriction"> /// The subgrid, on which the ReInit is performed /// </param> /// <param name="IncludingInterface"> /// !! Not yet functional !! /// True, if the subgrid contains the interface, this causes all external edges of the subgrid to be treated as boundaries /// False, for the rest of the domain, thus the flux to the adjacent cells wil be evaluated /// </param> /// <returns></returns> public void ReInitSingleStep(out double ChangeRate, SubGrid Restriction = null, bool IncludingInterface = true) { if (!IncludingInterface) { throw new NotImplementedException("Untested, not yet functional!"); } using (new FuncTrace()) { /// Init Residuals Residual.Clear(); Residual.Acc(1.0, Phi); OldPhi.Clear(); OldPhi.Acc(1.0, Phi); NewPhi.Clear(); NewPhi.Acc(1.0, Phi); CellMask RestrictionMask = Restriction == null ? null : Restriction.VolumeMask; //if (Control.Upwinding && UpdateDirection && IterationCounter % 10 == 0) { if (false && Control.Upwinding && UpdateDirection) { //if (Control.Upwinding && UpdateDirection) { UpdateBulkMatrix(Restriction); } UpdateDirection = false; // RHS part RHSField.CoordinateVector.Clear(); //Operator_RHS.Evaluate(NewPhi.Mapping, RHSField.Mapping); Operator_RHS.Evaluate(double.NaN, IncludingInterface ? Restriction : null, IncludingInterface ? SubGridBoundaryModes.BoundaryEdge : SubGridBoundaryModes.InnerEdge, ArrayTools.Cat(new DGField[] { Phi }, parameterFields, new DGField[] { RHSField })); #if DEBUG RHSField.CheckForNanOrInf(); #endif // solve // ===== double[] RHS = OpAffine.CloneAs(); RHS.ScaleV(-1.0); RHS.AccV(1.0, RHSField.CoordinateVector); SolverResult Result; if (Restriction != null) { SubRHS.Clear(); SubSolution.Clear(); SubRHS.AccV(1.0, RHS, default(int[]), SubVecIdx); SubSolution.AccV(1.0, NewPhi.CoordinateVector, default(int[]), SubVecIdx); Result = slv.Solve(SubSolution, SubRHS); NewPhi.Clear(RestrictionMask); NewPhi.CoordinateVector.AccV(1.0, SubSolution, SubVecIdx, default(int[])); } else { Result = slv.Solve(NewPhi.CoordinateVector, RHS); } #if Debug OpMatrix.SpMV(-1.0, NewPhi.CoordinateVector, 1.0, RHS); Console.WriteLine("LinearSolver: {0} Iterations, Converged={1}, Residual = {2} ", Result.NoOfIterations, Result.Converged, RHS.L2Norm()); #endif // Apply underrelaxation Phi.Clear(RestrictionMask); Phi.Acc(1 - underrelaxation, OldPhi, RestrictionMask); Phi.Acc(underrelaxation, NewPhi, RestrictionMask); Residual.Acc(-1.0, Phi, RestrictionMask); ChangeRate = Residual.L2Norm(RestrictionMask); //Calculate LevelSetGradient.Clear(); LevelSetGradient.Gradient(1.0, Phi, RestrictionMask); //LevelSetGradient.GradientByFlux(1.0, Phi); MeanLevelSetGradient.Clear(); MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient, RestrictionMask); if (Control.Upwinding) { //RestrictionMask.GetBitMask(); for (int i = 0; i < MeanLevelSetGradient.CoordinateVector.Length; i++) { NewDirection[i] = Math.Sign(MeanLevelSetGradient.CoordinateVector[i]); //NewDirection[i] = MeanLevelSetGradient.CoordinateVector[i]; OldDirection[i] -= NewDirection[i]; } double MaxDiff = OldDirection.L2Norm(); //if (MaxDiff > 1E-20 && IterationCounter % 10 == 0 ) { //if (MaxDiff > 1E-20) { // Console.WriteLine("Direction Values differ by {0}", MaxDiff); if (MaxDiff > 0.2) { //UpdateDirection = true; //Console.WriteLine("Direction Values differ by {0} => Updating ReInit-Matrix", MaxDiff); } ; //} //Console.WriteLine("HACK!!! Updating Upwind Matrix everytime!"); //UpdateDirection = true; // Reset Value OldDirection.Clear(); OldDirection.AccV(1.0, NewDirection); } } }
/// <summary> /// Solution of the system /// <see cref="LaplaceMtx"/>*<see cref="T"/> + <see cref="LaplaceAffine"/> = <see cref="RHS"/> /// using a black-box solver /// </summary> private void ClassicSolve(out double mintime, out double maxtime, out bool Converged, out int NoOfIter) { mintime = double.MaxValue; maxtime = double.MinValue; Converged = false; NoOfIter = int.MaxValue; for (int i = 0; i < base.Control.NoOfSolverRuns; i++) { // create sparse solver // -------------------- ISparseSolver ipSolver; LinearSolverCode solvercodes = LinearSolverCode.classic_pardiso; switch (solvercodes) { case LinearSolverCode.classic_pardiso: ipSolver = new ilPSP.LinSolvers.PARDISO.PARDISOSolver() { CacheFactorization = true, UseDoublePrecision = true }; break; case LinearSolverCode.classic_mumps: ipSolver = new ilPSP.LinSolvers.MUMPS.MUMPSSolver(); break; case LinearSolverCode.classic_cg: ipSolver = new ilPSP.LinSolvers.monkey.CG() { MaxIterations = 1000000, Tolerance = 1.0e-10, DevType = ilPSP.LinSolvers.monkey.DeviceType.Cuda }; break; default: throw new ArgumentException(); } ipSolver.DefineMatrix(LaplaceMtx); // call solver // ----------- T.Clear(); Console.WriteLine("RUN " + i + ": solving system..."); var RHSvec = RHS.CoordinateVector.ToArray(); BLAS.daxpy(RHSvec.Length, -1.0, this.LaplaceAffine, 1, RHSvec, 1); T.Clear(); SolverResult solRes = ipSolver.Solve(T.CoordinateVector, RHSvec); mintime = Math.Min(solRes.RunTime.TotalSeconds, mintime); maxtime = Math.Max(solRes.RunTime.TotalSeconds, maxtime); Converged = solRes.Converged; NoOfIter = solRes.NoOfIterations; Console.WriteLine("Pardiso phase 11: " + ilPSP.LinSolvers.PARDISO.PARDISOSolver.Phase_11.Elapsed.TotalSeconds); Console.WriteLine("Pardiso phase 22: " + ilPSP.LinSolvers.PARDISO.PARDISOSolver.Phase_22.Elapsed.TotalSeconds); Console.WriteLine("Pardiso phase 33: " + ilPSP.LinSolvers.PARDISO.PARDISOSolver.Phase_33.Elapsed.TotalSeconds); ipSolver.Dispose(); } }
/// <summary> /// Reinit on un-cut cells. /// </summary> /// <param name="Phi">The level set</param> /// <param name="ReInitSpecies">Cell mask wich is to be reinitialized</param> /// <param name="sign">Sign of the level set for this <paramref name="ReInitSpecies"/></param> /// <param name="_Accepted">CellMask which is taken as boundray values</param> /// <param name="GradPhi">LEvel Set gradient</param> /// <param name="callBack">A delegate, which might be called after the execution of the reinitialization</param> public void Reinitialize(SinglePhaseField Phi, CellMask ReInitSpecies, double sign, CellMask _Accepted, //ConventionalDGField[] ExtProperty, double[][] ExtPropertyMin, double[][] ExtPropertyMax, VectorField <SinglePhaseField> GradPhi, Action <int> callBack) // { using (new FuncTrace()) { Tracer.InstrumentationSwitch = false; // lots of tracing on calls acting on singe cells causes massive overhead (up to 5x slower). Stpw_total.Start(); SinglePhaseField DiffusionCoeff = new SinglePhaseField(new Basis(this.GridDat, 1), "DiffusionCoeff"); // check args and init // =================== /* * ExtVelSolver extVelSlv = null; * if(ExtProperty != null) { * if(ExtProperty.Length != ExtPropertyMin.Length) * throw new ArgumentException(); * if(ExtProperty.Length != ExtPropertyMax.Length) * throw new ArgumentException(); * * extVelSlv = new ExtVelSolver(ExtProperty[0].Basis); * } */ BitArray Acceped_Mutuable = _Accepted.GetBitMask().CloneAs(); BitArray Trial_Mutuable = ((_Accepted.AllNeighbourCells().Intersect(ReInitSpecies)).Except(_Accepted)).GetBitMask().CloneAs(); BitArray Recalc_Mutuable = Trial_Mutuable.CloneAs(); BitArray PosSpecies_Bitmask = ReInitSpecies.GetBitMask(); int J = this.GridDat.Cells.NoOfCells; int D = this.GridDat.SpatialDimension; int N = this.LevelSetBasis.Length; double _sign = sign >= 0 ? 1.0 : -1.0; double[] PhiAvg = m_PhiAvg; if (PhiAvg == null) { throw new ApplicationException(); } foreach (int jCell in _Accepted.ItemEnum) { PhiAvg[jCell] = Phi.GetMeanValue(jCell); } int NoOfNew; { var Neu = ReInitSpecies.Except(_Accepted); NoOfNew = Neu.NoOfItemsLocally; Phi.Clear(Neu); Phi.AccConstant(_sign, Neu); foreach (int jCell in Neu.ItemEnum) { PhiAvg[jCell] = 1.0e10; } } if (this.GridDat.MpiSize > 1) { throw new NotSupportedException("Currently not MPI parallel."); } for (int d = 0; d < this.GridDat.SpatialDimension; d++) { if (!GradPhi[d].Basis.Equals(Phi.Basis)) { throw new ArgumentException("Level-set and level-set gradient field should have the same DG basis."); // ein grad niedriger wrürde auch genügen... } } // perform marching... // =================== // update gradient for cut-cells GradPhi.Clear(_Accepted); GradPhi.Gradient(1.0, Phi, _Accepted); // marching loop../ int cnt = 0; while (true) { cnt++; CellMask Recalc = new CellMask(this.GridDat, Recalc_Mutuable); CellMask Accepted = new CellMask(this.GridDat, Acceped_Mutuable); CellMask Trial = new CellMask(this.GridDat, Trial_Mutuable); int NoOfTrial = Trial.NoOfItemsLocally; int NoOfAccpt = Accepted.NoOfItemsLocally; int NoOfRcalc = Recalc.NoOfItemsLocally; if (Trial.NoOfItemsLocally <= 0) { //Ploti(Recalc, Accepted, Trial, Phi, Phi_gradient, optEikonalOut, cnt); break; } // Local solver for all 'Recalc'-cells // -------------------------------------- if (Recalc.NoOfItemsLocally > 0) { this.LocalSolve(Accepted, Recalc, Phi, GradPhi, _sign, DiffusionCoeff); } // find the next cell to accept // ---------------------------- // get mean value in all cells foreach (int jCell in Recalc.ItemEnum) { PhiAvg[jCell] = Phi.GetMeanValue(jCell); Recalc_Mutuable[jCell] = false; } //Ploti(Recalc, Accepted, Trial, Phi, Phi_gradient, optEikonalOut, cnt); // find trial-cell with minimum average value // this should be done with heap-sort (see fast-marching algorithm) int jCellAccpt = int.MaxValue; double TrialMin = double.MaxValue; foreach (int jCell in Trial.ItemEnum) { if (PhiAvg[jCell] * _sign < TrialMin) { TrialMin = PhiAvg[jCell] * _sign; jCellAccpt = jCell; } } if (callBack != null) { callBack(cnt); } /* * // update the gradient * // ------------------- * * this.Stpw_gradientEval.Start(); * gradModule.GradientUpdate(jCellAccpt, Acceped_Mutuable, Phi, GradPhi); * this.Stpw_gradientEval.Stop(); * * /* * // solve for the extension properties * // ---------------------------------- * * if(ExtProperty != null) { * int[] Neight, dummy33; * GridDat.Cells.GetCellNeighbours(jCellAccpt, GridData.CellData.GetCellNeighbours_Mode.ViaEdges, out Neight, out dummy33); * * for(int iComp = 0; iComp < ExtProperty.Length; iComp++) { * * ExtPropertyMax[iComp][jCellAccpt] = -double.MaxValue; * ExtPropertyMin[iComp][jCellAccpt] = double.MaxValue; * * foreach(int jNeig in Neight) { * if(Acceped_Mutuable[jNeig]) { * ExtPropertyMax[iComp][jCellAccpt] = Math.Max(ExtPropertyMax[iComp][jCellAccpt], ExtPropertyMax[iComp][jNeig]); * ExtPropertyMin[iComp][jCellAccpt] = Math.Min(ExtPropertyMin[iComp][jCellAccpt], ExtPropertyMin[iComp][jNeig]); * } * } * * this.Stpw_extVelSolver.Start(); * extVelSlv.ExtVelSolve_Far(Phi, GradPhi, ExtProperty[iComp], ref ExtPropertyMin[iComp][jCellAccpt], ref ExtPropertyMax[iComp][jCellAccpt], jCellAccpt, Accepted, _sign); * this.Stpw_extVelSolver.Stop(); * } * } * /* * { * int[] Neight, dummy33; * GridDat.Cells.GetCellNeighbours(jCellAccpt, GridData.CellData.GetCellNeighbours_Mode.ViaEdges, out Neight, out dummy33); * foreach(int jNeig in Neight) { * if(Acceped_Mutuable[jNeig]) { * plotDependencyArrow(cnt, jCellAccpt, jNeig); * } * } * } */ // the mimium is moved to accepted // ------------------------------- Acceped_Mutuable[jCellAccpt] = true; Trial_Mutuable[jCellAccpt] = false; Recalc_Mutuable[jCellAccpt] = false; NoOfNew--; // recalc on all neighbours // ------------------------ int[] Neighs, dummy; this.GridDat.GetCellNeighbours(jCellAccpt, GetCellNeighbours_Mode.ViaEdges, out Neighs, out dummy); foreach (int jNeig in Neighs) { if (!Acceped_Mutuable[jNeig] && PosSpecies_Bitmask[jNeig]) { Trial_Mutuable[jNeig] = true; Recalc_Mutuable[jNeig] = true; } } } if (NoOfNew > 0) { throw new ArithmeticException("Unable to perform reinitialization for all requested cells - maybe they are not reachable from the initialy 'accepted' domain?"); } //PlottAlot("dependencies.csv"); Tracer.InstrumentationSwitch = true; Stpw_total.Stop(); } }
// Local Variables for Iteration // <summary> // Counter for Iteration Steps // </summary> //double OldResidual = double.MaxValue; //int divergencecounter = 0; ///// <summary> ///// Checks for Reaching Max. Number of Iterations and Divergence of Algorithm ///// </summary> ///// <param name="Residual">Change Rate of the Algorithm</param> ///// <returns>Reaching Max Iterations, Aborts when diverged</returns> //public bool CheckAbortCriteria(double Residual, int IterationCounter) { // if (Residual <= ConvergenceCriterion) { // Console.WriteLine("EllipticReInit converged after {0} Iterations ", IterationCounter); // return true; // } // if (Residual >= OldResidual) divergencecounter++; // else divergencecounter = 0; // if (IterationCounter >= MaxIteration) { // Console.WriteLine("Elliptic Reinit Max Iterations Reached"); // return true; // }; // if (divergencecounter > MaxIteration / 2) { // Console.WriteLine("Elliptic Reinit diverged - Aborting"); // throw new ApplicationException(); // } // OldResidual = Residual; // IterationCounter++; // return false; //} //bool PreviouslyOnSubgrid = false; /// <summary> /// Updates the Operator Matrix after level-set motion /// </summary> /// <param name="Restriction"> /// The subgrid, on which the ReInit is performed /// </param> /// <param name="IncludingInterface"> /// !! Not yet functional !! /// True, if the subgrid contains the interface, this causes all external edges of the subgrid to be treated as boundaries /// False, for the rest of the domain, thus the flux to the adjacent cells wil be evaluated /// </param> public void UpdateOperators(SubGrid Restriction = null, bool IncludingInterface = true) { if (!IncludingInterface) { throw new NotImplementedException("Untested, not yet functional!"); } using (new FuncTrace()) { //using (var slv = new ilPSP.LinSolvers.MUMPS.MUMPSSolver()) { //using (var slv = new ilPSP.LinSolvers.PARDISO.PARDISOSolver()) { //using (var slv = new ilPSP.LinSolvers.HYPRE.GMRES()) { if (Control.Upwinding) { OldPhi.Clear(); OldPhi.Acc(1.0, Phi); //Calculate LevelSetGradient.Clear(); LevelSetGradient.Gradient(1.0, Phi, Restriction?.VolumeMask); //LevelSetGradient.Gradient(1.0, Phi); //LevelSetGradient.GradientByFlux(1.0, Phi); MeanLevelSetGradient.Clear(); MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient, Restriction?.VolumeMask); //MeanLevelSetGradient.AccLaidBack(1.0, LevelSetGradient); } if (slv != null) { slv.Dispose(); } slv = Control.solverFactory(); OpMatrix_interface.Clear(); OpAffine_interface.Clear(); // Build the Quadrature-Scheme for the interface operator // Note: The HMF-Quadrature over a surface is formally a volume quadrature, since it uses the volume quadrature nodes. //XSpatialOperatorExtensions.ComputeMatrixEx(Operator_interface, ////Operator_interface.ComputeMatrixEx( // LevelSetTracker, // Phi.Mapping, // null, // Phi.Mapping, // OpMatrix_interface, // OpAffine_interface, // false, // 0, // false, // subGrid:Restriction, // whichSpc: LevelSetTracker.GetSpeciesId("A") // ); XSpatialOperatorMk2.XEvaluatorLinear mtxBuilder = Operator_interface.GetMatrixBuilder(LevelSetTracker, Phi.Mapping, null, Phi.Mapping); MultiphaseCellAgglomerator dummy = LevelSetTracker.GetAgglomerator(LevelSetTracker.SpeciesIdS.ToArray(), Phi.Basis.Degree * 2 + 2, 0.0); //mtxBuilder.SpeciesOperatorCoefficients[LevelSetTracker.GetSpeciesId("A")].CellLengthScales = dummy.CellLengthScales[LevelSetTracker.GetSpeciesId("A")]; mtxBuilder.CellLengthScales.Add(LevelSetTracker.GetSpeciesId("A"), dummy.CellLengthScales[LevelSetTracker.GetSpeciesId("A")]); mtxBuilder.time = 0; mtxBuilder.MPITtransceive = false; mtxBuilder.ComputeMatrix(OpMatrix_interface, OpAffine_interface); // Regenerate OpMatrix for subgrid -> adjacent cells must be trated as boundary if (Restriction != null) { OpMatrix_bulk.Clear(); OpAffine_bulk.Clear(); //Operator_bulk.ComputeMatrix( // Phi.Mapping, // parameterFields, // Phi.Mapping, // OpMatrix_bulk, OpAffine_bulk, // OnlyAffine: false, sgrd: Restriction); EdgeQuadratureScheme edgescheme; //if (Control.Upwinding) { // edgescheme = new EdgeQuadratureScheme(true, IncludingInterface ? Restriction.AllEdgesMask : null); //} //else { edgescheme = new EdgeQuadratureScheme(true, IncludingInterface ? Restriction.InnerEdgesMask : null); //} Operator_bulk.ComputeMatrixEx(Phi.Mapping, parameterFields, Phi.Mapping, OpMatrix_bulk, OpAffine_bulk, false, 0, edgeQuadScheme: edgescheme, volQuadScheme: new CellQuadratureScheme(true, IncludingInterface ? Restriction.VolumeMask : null) ); //PreviouslyOnSubgrid = true; } // recalculate full Matrix //else if (PreviouslyOnSubgrid) { else { OpMatrix_bulk.Clear(); OpAffine_bulk.Clear(); Operator_bulk.ComputeMatrixEx(Phi.Mapping, parameterFields, Phi.Mapping, OpMatrix_bulk, OpAffine_bulk, false, 0 ); //PreviouslyOnSubgrid = false; } /// Compose the Matrix /// This is symmetric due to the symmetry of the SIP and the penalty term OpMatrix.Clear(); OpMatrix.Acc(1.0, OpMatrix_bulk); OpMatrix.Acc(1.0, OpMatrix_interface); OpMatrix.AssumeSymmetric = !Control.Upwinding; //OpMatrix.AssumeSymmetric = false; /// Compose the RHS of the above operators. (-> Boundary Conditions) /// This does NOT include the Nonlinear RHS, which will be added later OpAffine.Clear(); OpAffine.AccV(1.0, OpAffine_bulk); OpAffine.AccV(1.0, OpAffine_interface); #if Debug ilPSP.Connectors.Matlab.BatchmodeConnector matlabConnector; matlabConnector = new BatchmodeConnector(); #endif if (Restriction != null) { SubVecIdx = Phi.Mapping.GetSubvectorIndices(Restriction, true, new int[] { 0 }); int L = SubVecIdx.Length; SubMatrix = new MsrMatrix(L); SubRHS = new double[L]; SubSolution = new double[L]; OpMatrix.AccSubMatrixTo(1.0, SubMatrix, SubVecIdx, default(int[]), SubVecIdx, default(int[])); slv.DefineMatrix(SubMatrix); #if Debug Console.WriteLine("ConditionNumber of ReInit-Matrix is " + SubMatrix.condest().ToString("E")); #endif } else { slv.DefineMatrix(OpMatrix); #if Debug Console.WriteLine("ConditionNumber of ReInit-Matrix is " + OpMatrix.condest().ToString("E")); #endif } } }