void NonLinearStatic :: updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d) // // updates some component, which is used by numerical method // to newly reached state. used mainly by numerical method // when new tangent stiffness is needed during finding // of new equilibrium stage. // { switch ( cmpn ) { case NonLinearLhs: if ( stiffMode == nls_tangentStiffness ) { stiffnessMatrix->zero(); // zero stiffness matrix #ifdef VERBOSE OOFEM_LOG_DEBUG("Assembling tangent stiffness matrix\n"); #endif this->assemble(*stiffnessMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), d); } else if ( ( stiffMode == nls_secantStiffness ) || ( stiffMode == nls_secantInitialStiffness && initFlag ) ) { #ifdef VERBOSE OOFEM_LOG_DEBUG("Assembling secant stiffness matrix\n"); #endif stiffnessMatrix->zero(); // zero stiffness matrix this->assemble(*stiffnessMatrix, tStep, TangentAssembler(SecantStiffness), EModelDefaultEquationNumbering(), d); initFlag = 0; } else if ( ( stiffMode == nls_elasticStiffness ) && ( initFlag || ( this->giveMetaStep( tStep->giveMetaStepNumber() )->giveFirstStepNumber() == tStep->giveNumber() ) ) ) { #ifdef VERBOSE OOFEM_LOG_DEBUG("Assembling elastic stiffness matrix\n"); #endif stiffnessMatrix->zero(); // zero stiffness matrix this->assemble(*stiffnessMatrix, tStep, TangentAssembler(ElasticStiffness), EModelDefaultEquationNumbering(), d); initFlag = 0; } else { // currently no action , this method is mainly intended to // assemble new tangent stiffness after each iteration // when secantStiffMode is on, we use the same stiffness // during iteration process } break; case InternalRhs: #ifdef VERBOSE OOFEM_LOG_DEBUG("Updating internal forces\n"); #endif // update internalForces and internalForcesEBENorm concurrently this->giveInternalForces(internalForces, true, d->giveNumber(), tStep); break; default: OOFEM_ERROR("Unknown Type of component."); } }
void PrescribedGradientBCPeriodic :: computeTangent(FloatMatrix &E, TimeStep *tStep) { EModelDefaultEquationNumbering fnum; DofIDEquationNumbering pnum(true, strain_id); EngngModel *rve = this->giveDomain()->giveEngngModel(); ///@todo Get this from engineering model std :: unique_ptr< SparseLinearSystemNM > solver( classFactory.createSparseLinSolver( ST_Petsc, this->domain, this->domain->giveEngngModel() ) ); // = rve->giveLinearSolver(); SparseMtrxType stype = solver->giveRecommendedMatrix(true); std :: unique_ptr< SparseMtrx > Kff( classFactory.createSparseMtrx( stype ) ); std :: unique_ptr< SparseMtrx > Kfp( classFactory.createSparseMtrx( stype ) ); std :: unique_ptr< SparseMtrx > Kpp( classFactory.createSparseMtrx( stype ) ); Kff->buildInternalStructure(rve, this->domain->giveNumber(), fnum); Kfp->buildInternalStructure(rve, this->domain->giveNumber(), fnum, pnum); Kpp->buildInternalStructure(rve, this->domain->giveNumber(), pnum); rve->assemble(*Kff, tStep, TangentAssembler(TangentStiffness), fnum, this->domain); rve->assemble(*Kfp, tStep, TangentAssembler(TangentStiffness), fnum, pnum, this->domain); rve->assemble(*Kpp, tStep, TangentAssembler(TangentStiffness), pnum, this->domain); int neq = Kfp->giveNumberOfRows(); int nsd = this->domain->giveNumberOfSpatialDimensions(); int ncomp = nsd * nsd; FloatMatrix grad_pert(ncomp, ncomp), rhs, sol(neq, ncomp); grad_pert.resize(ncomp, ncomp); // In fact, npeq should most likely equal ndev grad_pert.beUnitMatrix(); // Compute the solution to each of the pertubation of eps Kfp->times(grad_pert, rhs); solver->solve(*Kff, rhs, sol); // Compute the solution to each of the pertubation of eps FloatMatrix E_tmp; Kfp->timesT(sol, E_tmp); // Assuming symmetry of stiffness matrix // This is probably always zero, but for generality FloatMatrix tmpMat; Kpp->times(grad_pert, tmpMat); E_tmp.subtract(tmpMat); E_tmp.times( - 1.0 / ( this->domainSize(this->giveDomain(), this->set) + this->domainSize(this->giveDomain(), this->masterSet) )); E.resize(E_tmp.giveNumberOfRows(), E_tmp.giveNumberOfColumns()); if ( nsd == 3 ) { if ( E_tmp.giveNumberOfRows() == 6 ) { E.assemble(E_tmp, {1, 6, 5, 6, 2, 4, 5, 4, 3}); } else { E.assemble(E_tmp, {1, 9, 8, 6, 2, 7, 5, 4, 3}); } } else if ( nsd == 2 ) { E.assemble(E_tmp, {1, 4, 3, 2}); } else { E = E_tmp; } }
void TransientTransportProblem :: updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d) { // F(T) + C*dT/dt = Q // Linearized: // F(T^(k)) + K*a*dT_1 = Q - C * dT/dt^(k) - C/dt * dT_1 // Rearranged // (a*K + C/dt) * dT_1 = Q - (F(T^(k)) + C * dT/dt^(k)) // K_eff * dT_1 = Q - F_eff // Update: // T_1 += dT_1 ///@todo NRSolver should report when the solution changes instead of doing it this way. this->field->update(VM_Total, tStep, solution, EModelDefaultEquationNumbering()); ///@todo Need to reset the boundary conditions properly since some "update" is doing strange /// things such as applying the (wrong) boundary conditions. This call will be removed when that code can be removed. this->field->applyBoundaryCondition(tStep); if ( cmpn == InternalRhs ) { // F_eff = F(T^(k)) + C * dT/dt^(k) this->internalForces.zero(); this->assembleVector(this->internalForces, tStep, InternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1), & this->eNorm); this->updateSharedDofManagers(this->internalForces, EModelDefaultEquationNumbering(), InternalForcesExchangeTag); if ( lumped ) { // Note, inertia contribution cannot be computed on element level when lumped mass matrices are used. FloatArray oldSolution, vel; this->field->initialize(VM_Total, tStep->givePreviousStep(), oldSolution, EModelDefaultEquationNumbering()); vel.beDifferenceOf(solution, oldSolution); vel.times( 1./tStep->giveTimeIncrement() ); for ( int i = 0; i < vel.giveSize(); ++i ) { this->internalForces[i] += this->capacityDiag[i] * vel[i]; } } else { FloatArray tmp; this->assembleVector(this->internalForces, tStep, InertiaForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1), & tmp); this->eNorm.add(tmp); ///@todo Fix this, assembleVector shouldn't zero eNorm inside the functions. / Mikael } } else if ( cmpn == NonLinearLhs ) { // K_eff = (a*K + C/dt) if ( !this->keepTangent ) { this->effectiveMatrix->zero(); this->assemble( *effectiveMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), this->giveDomain(1) ); effectiveMatrix->times(alpha); if ( lumped ) { effectiveMatrix->addDiagonal(1./tStep->giveTimeIncrement(), capacityDiag); } else { effectiveMatrix->add(1./tStep->giveTimeIncrement(), *capacityMatrix); } } } else { OOFEM_ERROR("Unknown component"); } }
void StaticStructural :: updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d) { if ( cmpn == InternalRhs ) { // Updates the solution in case it has changed ///@todo NRSolver should report when the solution changes instead of doing it this way. this->field->update(VM_Total, tStep, this->solution, EModelDefaultEquationNumbering()); this->field->applyBoundaryCondition(tStep);///@todo Temporary hack to override the incorrect vavues that is set by "update" above. Remove this when that is fixed. this->internalForces.zero(); this->assembleVector(this->internalForces, tStep, InternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), d, & this->eNorm); this->updateSharedDofManagers(this->internalForces, EModelDefaultEquationNumbering(), InternalForcesExchangeTag); internalVarUpdateStamp = tStep->giveSolutionStateCounter(); // Hack for linearstatic } else if ( cmpn == NonLinearLhs ) { this->stiffnessMatrix->zero(); this->assemble(*this->stiffnessMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), d); } else { OOFEM_ERROR("Unknown component"); } }
void StokesFlow :: updateComponent(TimeStep *tStep, NumericalCmpn cmpn, Domain *d) { velocityPressureField->update(VM_Total, tStep, solutionVector, EModelDefaultEquationNumbering()); // update element stabilization for ( auto &elem : d->giveElements() ) { static_cast< FMElement * >( elem.get() )->updateStabilizationCoeffs(tStep); } if ( cmpn == InternalRhs ) { this->internalForces.zero(); this->assembleVector(this->internalForces, tStep, InternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), d, & this->eNorm); this->updateSharedDofManagers(this->internalForces, EModelDefaultEquationNumbering(), InternalForcesExchangeTag); return; } else if ( cmpn == NonLinearLhs ) { this->stiffnessMatrix->zero(); this->assemble(*stiffnessMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), d); return; } else { OOFEM_ERROR("Unknown component"); } }
void FreeWarping :: solveYourselfAt(TimeStep *tStep) { // // creates system of governing eq's and solves them at given time step // // first assemble problem at current time step if ( initFlag ) { #ifdef VERBOSE OOFEM_LOG_DEBUG("Assembling stiffness matrix\n"); #endif // // first step assemble stiffness Matrix // stiffnessMatrix = classFactory.createSparseMtrx(sparseMtrxType); if ( stiffnessMatrix == NULL ) { OOFEM_ERROR("sparse matrix creation failed"); } stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() ); this->assemble( *stiffnessMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), this->giveDomain(1) ); initFlag = 0; } #ifdef VERBOSE OOFEM_LOG_DEBUG("Assembling load\n"); #endif // // allocate space for displacementVector // displacementVector.resize( this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ) ); displacementVector.zero(); // // assembling the load vector // loadVector.resize( this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ) ); loadVector.zero(); this->assembleVector( loadVector, tStep, ExternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1) ); // // internal forces (from Dirichlet b.c's, or thermal expansion, etc.) // // no such forces exist for the free warping problem /* FloatArray internalForces( this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ) ); internalForces.zero(); this->assembleVector( internalForces, tStep, InternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1) ); loadVector.subtract(internalForces); */ this->updateSharedDofManagers(loadVector, EModelDefaultEquationNumbering(), ReactionExchangeTag); // // set-up numerical model // this->giveNumericalMethod( this->giveMetaStep( tStep->giveMetaStepNumber() ) ); // // call numerical model to solve arose problem // #ifdef VERBOSE OOFEM_LOG_INFO("\n\nSolving ...\n\n"); #endif NM_Status s = nMethod->solve(*stiffnessMatrix, loadVector, displacementVector); if ( !( s & NM_Success ) ) { OOFEM_ERROR("No success in solving system."); } tStep->incrementStateCounter(); // update solution state counter }
void TransportGradientPeriodic :: computeTangent(FloatMatrix &k, TimeStep *tStep) { EModelDefaultEquationNumbering fnum; //DofIDEquationNumbering pnum(true, this->grad_ids); EModelDefaultPrescribedEquationNumbering pnum; int nsd = this->domain->giveNumberOfSpatialDimensions(); EngngModel *rve = this->giveDomain()->giveEngngModel(); ///@todo Get this from engineering model std :: unique_ptr< SparseLinearSystemNM > solver( classFactory.createSparseLinSolver( ST_Petsc, this->domain, this->domain->giveEngngModel() ) ); // = rve->giveLinearSolver(); SparseMtrxType stype = solver->giveRecommendedMatrix(true); std :: unique_ptr< SparseMtrx > Kff( classFactory.createSparseMtrx( stype ) ); std :: unique_ptr< SparseMtrx > Kfp( classFactory.createSparseMtrx( stype ) ); std :: unique_ptr< SparseMtrx > Kpp( classFactory.createSparseMtrx( stype ) ); Kff->buildInternalStructure(rve, this->domain->giveNumber(), fnum); int neq = Kff->giveNumberOfRows(); Kfp->buildInternalStructure(rve, this->domain->giveNumber(), fnum, pnum); Kpp->buildInternalStructure(rve, this->domain->giveNumber(), pnum); //Kfp->buildInternalStructure(rve, neq, nsd, {}, {}); //Kpp->buildInternalStructure(rve, nsd, nsd, {}, {}); #if 0 rve->assemble(*Kff, tStep, TangentAssembler(TangentStiffness), fnum, this->domain); rve->assemble(*Kfp, tStep, TangentAssembler(TangentStiffness), fnum, pnum, this->domain); rve->assemble(*Kpp, tStep, TangentAssembler(TangentStiffness), pnum, this->domain); #else auto ma = TangentAssembler(TangentStiffness); IntArray floc, ploc; FloatMatrix mat, R; int nelem = domain->giveNumberOfElements(); #ifdef _OPENMP #pragma omp parallel for shared(Kff, Kfp, Kpp) private(mat, R, floc, ploc) #endif for ( int ielem = 1; ielem <= nelem; ielem++ ) { Element *element = domain->giveElement(ielem); // skip remote elements (these are used as mirrors of remote elements on other domains // when nonlocal constitutive models are used. They introduction is necessary to // allow local averaging on domains without fine grain communication between domains). if ( element->giveParallelMode() == Element_remote || !element->isActivated(tStep) ) { continue; } ma.matrixFromElement(mat, *element, tStep); if ( mat.isNotEmpty() ) { ma.locationFromElement(floc, *element, fnum); ma.locationFromElement(ploc, *element, pnum); ///@todo This rotation matrix is not flexible enough.. it can only work with full size matrices and doesn't allow for flexibility in the matrixassembler. if ( element->giveRotationMatrix(R) ) { mat.rotatedWith(R); } #ifdef _OPENMP #pragma omp critical #endif { Kff->assemble(floc, mat); Kfp->assemble(floc, ploc, mat); Kpp->assemble(ploc, mat); } } } Kff->assembleBegin(); Kfp->assembleBegin(); Kpp->assembleBegin(); Kff->assembleEnd(); Kfp->assembleEnd(); Kpp->assembleEnd(); #endif FloatMatrix grad_pert(nsd, nsd), rhs, sol(neq, nsd); grad_pert.resize(nsd, nsd); grad_pert.beUnitMatrix(); // Workaround since the matrix size is inflexible with custom dof numbering (so far, planned to be fixed). IntArray grad_loc; this->grad->giveLocationArray(this->grad_ids, grad_loc, pnum); FloatMatrix pert(Kpp->giveNumberOfRows(), nsd); pert.assemble(grad_pert, grad_loc, {1,2,3}); //pert.printYourself("pert"); //printf("Kfp = %d x %d\n", Kfp->giveNumberOfRows(), Kfp->giveNumberOfColumns()); //printf("Kff = %d x %d\n", Kff->giveNumberOfRows(), Kff->giveNumberOfColumns()); //printf("Kpp = %d x %d\n", Kpp->giveNumberOfRows(), Kpp->giveNumberOfColumns()); // Compute the solution to each of the pertubation of eps Kfp->times(pert, rhs); //rhs.printYourself("rhs"); // Initial guess (Taylor assumption) helps KSP-iterations for ( auto &n : domain->giveDofManagers() ) { int k1 = n->giveDofWithID( this->dofs(0) )->__giveEquationNumber(); if ( k1 ) { FloatArray *coords = n->giveCoordinates(); for ( int i = 1; i <= nsd; ++i ) { sol.at(k1, i) = -(coords->at(i) - mCenterCoord.at(i)); } } } if ( solver->solve(*Kff, rhs, sol) & NM_NoSuccess ) { OOFEM_ERROR("Failed to solve Kff"); } // Compute the solution to each of the pertubation of eps Kfp->timesT(sol, k); // Assuming symmetry of stiffness matrix // This is probably always zero, but for generality FloatMatrix tmpMat; Kpp->times(pert, tmpMat); k.subtract(tmpMat); k.times( - 1.0 / ( this->domainSize(this->giveDomain(), this->set) + this->domainSize(this->giveDomain(), this->masterSet) )); // Temp workaround on sizing issue mentioned above: FloatMatrix k2 = k; k.beSubMatrixOf(k2, grad_loc, {1,2,3}); }
void IncrementalLinearStatic :: solveYourselfAt(TimeStep *tStep) { Domain *d = this->giveDomain(1); // Creates system of governing eq's and solves them at given time step // >>> beginning PH // The following piece of code updates assignment of boundary conditions to dofs // (this allows to have multiple boundary conditions assigned to one dof // which can be arbitrarily turned on and off in time) // Almost the entire section has been copied from domain.C std :: vector< std :: map< int, int > > dof_bc( d->giveNumberOfDofManagers() ); for ( int i = 1; i <= d->giveNumberOfBoundaryConditions(); ++i ) { GeneralBoundaryCondition *gbc = d->giveBc(i); if ( gbc->isImposed(tStep) ) { if ( gbc->giveSetNumber() > 0 ) { ///@todo This will eventually not be optional. // Loop over nodes in set and store the bc number in each dof. Set *set = d->giveSet( gbc->giveSetNumber() ); ActiveBoundaryCondition *active_bc = dynamic_cast< ActiveBoundaryCondition * >(gbc); BoundaryCondition *bc = dynamic_cast< BoundaryCondition * >(gbc); if ( bc || ( active_bc && active_bc->requiresActiveDofs() ) ) { const IntArray &appliedDofs = gbc->giveDofIDs(); const IntArray &nodes = set->giveNodeList(); for ( int inode = 1; inode <= nodes.giveSize(); ++inode ) { for ( int idof = 1; idof <= appliedDofs.giveSize(); ++idof ) { if ( dof_bc [ nodes.at(inode) - 1 ].find( appliedDofs.at(idof) ) == dof_bc [ nodes.at(inode) - 1 ].end() ) { // is empty dof_bc [ nodes.at(inode) - 1 ] [ appliedDofs.at(idof) ] = i; DofManager * dofman = d->giveDofManager( nodes.at(inode) ); Dof * dof = dofman->giveDofWithID( appliedDofs.at(idof) ); dof->setBcId(i); } else { // another bc has been already prescribed at this time step to this dof OOFEM_WARNING("More than one boundary condition assigned at time %f to node %d dof %d. Considering boundary condition %d", tStep->giveTargetTime(), nodes.at(inode), appliedDofs.at(idof), dof_bc [ nodes.at(inode) - 1 ] [appliedDofs.at(idof)] ); } } } } } } } // to get proper number of equations this->forceEquationNumbering(); // <<< end PH // Initiates the total displacement to zero. if ( tStep->isTheFirstStep() ) { for ( auto &dofman : d->giveDofManagers() ) { for ( Dof *dof: *dofman ) { dof->updateUnknownsDictionary(tStep->givePreviousStep(), VM_Total, 0.); dof->updateUnknownsDictionary(tStep, VM_Total, 0.); } } for ( auto &bc : d->giveBcs() ) { ActiveBoundaryCondition *abc; if ( ( abc = dynamic_cast< ActiveBoundaryCondition * >(bc.get()) ) ) { int ndman = abc->giveNumberOfInternalDofManagers(); for ( int i = 1; i <= ndman; i++ ) { DofManager *dofman = abc->giveInternalDofManager(i); for ( Dof *dof: *dofman ) { dof->updateUnknownsDictionary(tStep->givePreviousStep(), VM_Total, 0.); dof->updateUnknownsDictionary(tStep, VM_Total, 0.); } } } } } // Apply dirichlet b.c's on total values for ( auto &dofman : d->giveDofManagers() ) { for ( Dof *dof: *dofman ) { double tot = dof->giveUnknown( VM_Total, tStep->givePreviousStep() ); if ( dof->hasBc(tStep) ) { tot += dof->giveBcValue(VM_Incremental, tStep); } dof->updateUnknownsDictionary(tStep, VM_Total, tot); } } int neq = this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ); #ifdef VERBOSE OOFEM_LOG_RELEVANT("Solving [step number %8d, time %15e, equations %d]\n", tStep->giveNumber(), tStep->giveTargetTime(), neq); #endif if ( neq == 0 ) { // Allows for fully prescribed/empty problems. return; } incrementOfDisplacementVector.resize(neq); incrementOfDisplacementVector.zero(); #ifdef VERBOSE OOFEM_LOG_INFO("Assembling load\n"); #endif // Assembling the element part of load vector internalLoadVector.resize(neq); internalLoadVector.zero(); this->assembleVector( internalLoadVector, tStep, InternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1) ); loadVector.resize(neq); loadVector.zero(); this->assembleVector( loadVector, tStep, ExternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1) ); loadVector.subtract(internalLoadVector); this->updateSharedDofManagers(loadVector, EModelDefaultEquationNumbering(), ReactionExchangeTag); #ifdef VERBOSE OOFEM_LOG_INFO("Assembling stiffness matrix\n"); #endif stiffnessMatrix.reset( classFactory.createSparseMtrx(sparseMtrxType) ); if ( !stiffnessMatrix ) { OOFEM_ERROR("sparse matrix creation failed"); } stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() ); stiffnessMatrix->zero(); this->assemble( *stiffnessMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), this->giveDomain(1) ); #ifdef VERBOSE OOFEM_LOG_INFO("Solving ...\n"); #endif this->giveNumericalMethod( this->giveCurrentMetaStep() ); NM_Status s = nMethod->solve(*stiffnessMatrix, loadVector, incrementOfDisplacementVector); if ( !( s & NM_Success ) ) { OOFEM_ERROR("No success in solving system."); } }
void TransientTransportProblem :: solveYourselfAt(TimeStep *tStep) { Domain *d = this->giveDomain(1); int neq = this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ); field->advanceSolution(tStep); field->applyBoundaryCondition(tStep); if ( tStep->isTheFirstStep() ) { this->applyIC(); } field->initialize(VM_Total, tStep, solution, EModelDefaultEquationNumbering()); if ( !effectiveMatrix ) { effectiveMatrix.reset( classFactory.createSparseMtrx(sparseMtrxType) ); effectiveMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() ); if ( lumped ) { capacityDiag.resize(neq); this->assembleVector( capacityDiag, tStep, LumpedMassVectorAssembler(), VM_Total, EModelDefaultEquationNumbering(), d ); } else { capacityMatrix.reset( classFactory.createSparseMtrx(sparseMtrxType) ); capacityMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() ); this->assemble( *capacityMatrix, tStep, MassMatrixAssembler(), EModelDefaultEquationNumbering(), d ); } if ( this->keepTangent ) { this->assemble( *effectiveMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), d ); effectiveMatrix->times(alpha); if ( lumped ) { effectiveMatrix->addDiagonal(1./tStep->giveTimeIncrement(), capacityDiag); } else { effectiveMatrix->add(1./tStep->giveTimeIncrement(), *capacityMatrix); } } } OOFEM_LOG_INFO("Assembling external forces\n"); FloatArray externalForces(neq); externalForces.zero(); this->assembleVector( externalForces, tStep, ExternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), d ); this->updateSharedDofManagers(externalForces, EModelDefaultEquationNumbering(), LoadExchangeTag); // set-up numerical method this->giveNumericalMethod( this->giveCurrentMetaStep() ); OOFEM_LOG_INFO("Solving for %d unknowns...\n", neq); internalForces.resize(neq); FloatArray incrementOfSolution; double loadLevel; int currentIterations; this->nMethod->solve(*this->effectiveMatrix, externalForces, NULL, // ignore this->solution, incrementOfSolution, this->internalForces, this->eNorm, loadLevel, // ignore SparseNonLinearSystemNM :: rlm_total, // ignore currentIterations, // ignore tStep); }
void DDLinearStatic :: solveYourselfAt(TimeStep *tStep) { /** * Perform DD * Currently this is completely specific, but needs to be generalized * by perhaps creating an DD_interface class to do set of functions */ /**************************************************************************/ int NumberOfDomains = this->giveNumberOfDomains(); // Loop through all the giveNumberOfDomains for (int i = 1; i<= NumberOfDomains; i++) { //Domain *domain = this->giveDomain(i); /// Create a spatial localizer which in effect has services for locating points in element etc. //SpatialLocalizer *sl = domain->giveSpatialLocalizer(); // Perform DD solution at this time step, solution will depend on quantities in the input file /// @todo these have to be defined and initialized earlier /// Each domain for the DD case can have only one material... throw error otherwise //// Also the DD_domains should be intialized with these materials properties during input //// Here i am just using values from input file for algorithmic convenience /* dd::OofemInterface * interface = new dd::OofemInterface(this); */ dd::Domain dd_domain(70e-3, 0.3, NULL); dd::SlipSystem ss0 = dd::SlipSystem(0.0, 0.25e-3); dd_domain.addSlipSystem(&ss0); dd::SlipPlane sp0 = dd::SlipPlane(&dd_domain, &ss0, 0.0); dd::ObstaclePoint o0 = dd::ObstaclePoint(&dd_domain, &sp0, -0.25, 20.0e3); dd::ObstaclePoint o1 = dd::ObstaclePoint(&dd_domain, &sp0, 0.25, 20.0e3); double e = dd_domain.getModulus(); double nu = dd_domain.getPassionsRatio(); double mu = e / (2. * ( 1. + nu)); double fact = mu * ss0.getBurgersMagnitude() / ( 2 * M_PI * (1. - nu)); dd::SourcePoint s1 = dd::SourcePoint(&dd_domain, &sp0, 0, 25e-6, fact / 25e-6); for( dd_domain.dtNo = 1; dd_domain.dtNo < dd_domain.dtNomax; dd_domain.dtNo++) { std::cerr << "Total dislocs in domain: " << sp0.getContainer<dd::DislocationPoint>().size() << "\n"; std::cerr << "Dislocs: " << sp0.dumpToString<dd::DislocationPoint>() << "\n"; std::cerr.flush(); dd_domain.updateForceCaches(); for(auto point : dd_domain.getContainer<dd::DislocationPoint>()) { dd::Vector<2> force, forceGradient; dd::Vector<3> stress; force = dd::Vector<2>({0.0,0.0}); //point->sumCaches(force, forceGradient, stress); force = point->cachedForce(); stress = point->cachedStress(); std::cout << "Cached Force at " << point->slipPlanePosition() << ": " << point->getBurgersSign() << " " << force[0] << " " << stress[2] << " " << point->slipPlanePosition() << "\n"; } dd_domain.updateForceCaches(); dd_domain.moveDislocations(1.0e-11, 1.0e-18); s1.spawn(1, 5); /* for(int bcNo = 1; bcNo <= giveDomain(i)->giveNumberOfBoundaryConditions(); bcNo++) { ManualBoundaryCondition * bc = dynamic_cast<ManualBoundaryCondition *>(giveDomain(i)->giveBc(bcNo)); if(bc == nullptr || bc->giveType() != DirichletBT) { continue; } dd::Vector<2> bcContribution; Domain * d = bc->giveDomain(); Set * set = d->giveSet(bc->giveSetNumber()); for(int nodeNo : set->giveNodeList()) { Node * node = static_cast<Node *>(d->giveDofManager(nodeNo)); for (auto &dofid : bc->giveDofIDs()) { Dof * dof = node->giveDofWithID(dofid); interface->giveNodalBcContribution(node, bcContribution); // TODO: Determine the dimensions without pointer checking double toAdd; if(dof->giveDofID() == D_u) { toAdd = bcContribution[0]; } else if(dof->giveDofID() == D_v) { toAdd = bcContribution[1]; } else { OOFEM_ERROR("DOF must be x-disp or y-disp"); } bc->addManualValue(dof, toAdd); } } std::cout << "BC Contribution: " << bcContribution[0] << " " << bcContribution[1] << "\n"; } */ //delete interface; } // end dtNo loop } /**************************************************************************/ // // creates system of governing eq's and solves them at given time step // // first assemble problem at current time step if ( initFlag ) { #ifdef VERBOSE OOFEM_LOG_DEBUG("Assembling stiffness matrix\n"); #endif // // first step assemble stiffness Matrix // stiffnessMatrix.reset( classFactory.createSparseMtrx(sparseMtrxType) ); if ( !stiffnessMatrix ) { OOFEM_ERROR("sparse matrix creation failed"); } stiffnessMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() ); this->assemble( *stiffnessMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), this->giveDomain(1) ); initFlag = 0; } #ifdef VERBOSE OOFEM_LOG_DEBUG("Assembling load\n"); #endif // // allocate space for displacementVector // displacementVector.resize( this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ) ); displacementVector.zero(); // // assembling the load vector // loadVector.resize( this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ) ); loadVector.zero(); this->assembleVector( loadVector, tStep, ExternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1) ); // // internal forces (from Dirichlet b.c's, or thermal expansion, etc.) // FloatArray internalForces( this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ) ); internalForces.zero(); this->assembleVector( internalForces, tStep, InternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), this->giveDomain(1) ); loadVector.subtract(internalForces); this->updateSharedDofManagers(loadVector, EModelDefaultEquationNumbering(), ReactionExchangeTag); // // set-up numerical model // this->giveNumericalMethod( this->giveMetaStep( tStep->giveMetaStepNumber() ) ); // // call numerical model to solve arose problem // #ifdef VERBOSE OOFEM_LOG_INFO("\n\nSolving ...\n\n"); #endif NM_Status s = nMethod->solve(*stiffnessMatrix, loadVector, displacementVector); if ( !( s & NM_Success ) ) { OOFEM_ERROR("No success in solving system."); } tStep->incrementStateCounter(); // update solution state counter }
void TransientTransportProblem :: solveYourselfAt(TimeStep *tStep) { Domain *d = this->giveDomain(1); int neq = this->giveNumberOfDomainEquations( 1, EModelDefaultEquationNumbering() ); if ( tStep->isTheFirstStep() ) { this->applyIC(); } field->advanceSolution(tStep); #if 1 // This is what advanceSolution should be doing, but it can't be there yet // (backwards compatibility issues due to inconsistencies in other solvers). TimeStep *prev = tStep->givePreviousStep(); for ( auto &dman : d->giveDofManagers() ) { static_cast< DofDistributedPrimaryField* >(field.get())->setInitialGuess(*dman, tStep, prev); } for ( auto &elem : d->giveElements() ) { int ndman = elem->giveNumberOfInternalDofManagers(); for ( int i = 1; i <= ndman; i++ ) { static_cast< DofDistributedPrimaryField* >(field.get())->setInitialGuess(*elem->giveInternalDofManager(i), tStep, prev); } } for ( auto &bc : d->giveBcs() ) { int ndman = bc->giveNumberOfInternalDofManagers(); for ( int i = 1; i <= ndman; i++ ) { static_cast< DofDistributedPrimaryField* >(field.get())->setInitialGuess(*bc->giveInternalDofManager(i), tStep, prev); } } #endif field->applyBoundaryCondition(tStep); field->initialize(VM_Total, tStep, solution, EModelDefaultEquationNumbering()); if ( !effectiveMatrix ) { effectiveMatrix.reset( classFactory.createSparseMtrx(sparseMtrxType) ); effectiveMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() ); if ( lumped ) { capacityDiag.resize(neq); this->assembleVector( capacityDiag, tStep, LumpedMassVectorAssembler(), VM_Total, EModelDefaultEquationNumbering(), d ); } else { capacityMatrix.reset( classFactory.createSparseMtrx(sparseMtrxType) ); capacityMatrix->buildInternalStructure( this, 1, EModelDefaultEquationNumbering() ); this->assemble( *capacityMatrix, tStep, MassMatrixAssembler(), EModelDefaultEquationNumbering(), d ); } if ( this->keepTangent ) { this->assemble( *effectiveMatrix, tStep, TangentAssembler(TangentStiffness), EModelDefaultEquationNumbering(), d ); effectiveMatrix->times(alpha); if ( lumped ) { effectiveMatrix->addDiagonal(1./tStep->giveTimeIncrement(), capacityDiag); } else { effectiveMatrix->add(1./tStep->giveTimeIncrement(), *capacityMatrix); } } } OOFEM_LOG_INFO("Assembling external forces\n"); FloatArray externalForces(neq); externalForces.zero(); this->assembleVector( externalForces, tStep, ExternalForceAssembler(), VM_Total, EModelDefaultEquationNumbering(), d ); this->updateSharedDofManagers(externalForces, EModelDefaultEquationNumbering(), LoadExchangeTag); // set-up numerical method this->giveNumericalMethod( this->giveCurrentMetaStep() ); OOFEM_LOG_INFO("Solving for %d unknowns...\n", neq); internalForces.resize(neq); FloatArray incrementOfSolution; double loadLevel; int currentIterations; this->nMethod->solve(*this->effectiveMatrix, externalForces, NULL, // ignore NULL, this->solution, incrementOfSolution, this->internalForces, this->eNorm, loadLevel, // ignore SparseNonLinearSystemNM :: rlm_total, // ignore currentIterations, // ignore tStep); }