Changeset 18350
- Timestamp:
- 08/08/14 13:52:03 (11 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/solutionsequences/solutionsequence_fct.cpp
r18349 r18350 112 112 *pLHS = LHS; 113 113 }/*}}}*/ 114 void CreateRHS(Vec* pRHS,Mat K,Mat D,Vec Ml,Vec u,IssmDouble theta,IssmDouble deltat,IssmDouble dmax,FemModel* femmodel,int configuration_type){/*{{{*/ 115 /*Create Left Hand side of Lower order solution 116 * 117 * RHS = [ML + (1 − theta) deltaT L^n] u^n 118 * 119 * where L = K + D 120 * 121 */ 122 123 /*Intermediaries*/ 124 Vec Ku = NULL; 125 Vec Du = NULL; 126 Vec RHS = NULL; 127 int dof; 128 IssmDouble d; 129 130 /*Initialize vectors*/ 131 VecDuplicate(u,&Ku); 132 VecDuplicate(u,&Du); 133 VecDuplicate(u,&RHS); 134 135 /*Create RHS = M*u + (1-theta)*deltat*K*u + (1-theta)*deltat*D*u*/ 136 MatMult(K,u,Ku); 137 MatMult(D,u,Du); 138 VecPointwiseMult(RHS,Ml,u); 139 VecAXPBYPCZ(RHS,(1-theta)*deltat,(1-theta)*deltat,1,Ku,Du); 140 VecFree(&Ku); 141 VecFree(&Du); 142 143 /*Penalize Dirichlet boundary*/ 144 for(int i=0;i<femmodel->constraints->Size();i++){ 145 Constraint* constraint=(Constraint*)femmodel->constraints->GetObjectByOffset(i); 146 if(constraint->InAnalysis(configuration_type)){ 147 constraint->PenaltyDofAndValue(&dof,&d,femmodel->nodes,femmodel->parameters); 148 d = d*dmax; 149 if(dof!=-1){ 150 VecSetValues(RHS,1,&dof,(const double*)&d,INSERT_VALUES); 151 } 152 } 153 } 154 VecAssemblyBegin(RHS); 155 VecAssemblyEnd( RHS); 156 157 /*Assign output pointer*/ 158 *pRHS = RHS; 159 160 }/*}}}*/ 114 161 #endif 115 162 void solutionsequence_fct(FemModel* femmodel){ … … 121 168 Vector<IssmDouble>* ug = NULL; 122 169 Vector<IssmDouble>* uf = NULL; 123 Vector<IssmDouble>* ys = NULL; 124 125 IssmDouble theta,deltat,dmax=0.; 126 int dof,ncols,ncols2,ncols3,rstart,rend; 127 int configuration_type,analysis_type; 128 double d,mi; 170 171 IssmDouble theta,deltat,dmax; 172 int dof,ncols,ncols2,rstart,rend; 173 int configuration_type; 174 double d; 129 175 int* cols = NULL; 130 176 int* cols2 = NULL; 131 int* cols3 = NULL;132 177 double* vals = NULL; 133 178 double* vals2 = NULL; 134 double* vals3 = NULL;135 179 136 180 /*Create analysis*/ … … 140 184 femmodel->parameters->FindParam(&deltat,TimesteppingTimeStepEnum); 141 185 femmodel->parameters->FindParam(&configuration_type,ConfigurationTypeEnum); 142 femmodel->parameters->FindParam(&analysis_type,AnalysisTypeEnum);143 186 femmodel->UpdateConstraintsx(); 144 187 theta = 0.5; … … 152 195 Mat D_petsc = NULL; 153 196 Mat LHS = NULL; 197 Vec RHS = NULL; 198 Vec u = NULL; 154 199 Mat K_petsc = K->pmatrix->matrix; 155 200 Vec Ml_petsc = Ml->pvector->vector; … … 167 212 Reducevectorgtofx(&uf, ug, femmodel->nodes,femmodel->parameters); 168 213 delete ug; 169 Vec u = uf->pvector->vector; 170 Vec Ku = NULL; 171 Vec Du = NULL; 172 Vec RHS = NULL; 173 VecDuplicate(u,&Ku); 174 VecDuplicate(u,&Du); 175 VecDuplicate(u,&RHS); 176 MatMult(K_petsc,u,Ku); 177 MatMult(D_petsc,u,Du); 178 VecPointwiseMult(RHS,Ml_petsc,u); 179 VecAXPBYPCZ(RHS,(1-theta)*deltat,(1-theta)*deltat,1,Ku,Du);// RHS = M*u + (1-theta)*deltat*K*u + (1-theta)*deltat*D*u 180 VecFree(&Ku); 181 VecFree(&Du); 214 CreateRHS(&RHS,K_petsc,D_petsc,Ml_petsc,uf->pvector->vector,theta,deltat,dmax,femmodel,configuration_type); 182 215 delete uf; 183 184 /*Penalize Dirichlet boundary*/185 for(int i=0;i<femmodel->constraints->Size();i++){186 Constraint* constraint=(Constraint*)femmodel->constraints->GetObjectByOffset(i);187 if(constraint->InAnalysis(analysis_type)){188 constraint->PenaltyDofAndValue(&dof,&d,femmodel->nodes,femmodel->parameters);189 d = d*dmax;190 if(dof!=-1){191 VecSetValues(RHS,1,&dof,(const double*)&d,INSERT_VALUES);192 }193 }194 }195 VecAssemblyBegin(RHS);196 VecAssemblyEnd( RHS);197 216 198 217 /*Go solve!*/ … … 325 344 VecFree(&u); 326 345 327 Mergesolutionfromftogx(&ug, uf,ys,femmodel->nodes,femmodel->parameters);delete uf; delete ys; 328 InputUpdateFromSolutionx(femmodel,ug); 329 delete ug; 346 InputUpdateFromSolutionx(femmodel,uf); 347 delete uf; 330 348 331 349 #else
Note:
See TracChangeset
for help on using the changeset viewer.