Changeset 4013
- Timestamp:
- 06/03/10 12:09:20 (15 years ago)
- Location:
- issm/trunk/src/c
- Files:
-
- 1 added
- 14 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk/src/c/EnumDefinitions/EnumDefinitions.h
r3982 r4013 46 46 //slope 47 47 SurfaceAnalysisEnum, 48 Slope computeAnalysisEnum,48 SlopeComputeAnalysisEnum, 49 49 BedXAnalysisEnum, 50 50 BedYAnalysisEnum, … … 154 154 ArtDiffEnum, 155 155 BedEnum, 156 BedSlope xEnum,157 BedSlope yEnum,156 BedSlopeXEnum, 157 BedSlopeYEnum, 158 158 BoundaryEnum, 159 159 CmMaxDmpSlopeEnum, … … 204 204 SurfaceAreaEnum, 205 205 SurfaceEnum, 206 SurfaceSlope xEnum,207 SurfaceSlope yEnum,206 SurfaceSlopeXEnum, 207 SurfaceSlopeYEnum, 208 208 TemperatureEnum, 209 209 TemperatureAverageEnum, -
issm/trunk/src/c/Makefile.am
r4002 r4013 951 951 ./solutions/transient_core_3d.cpp\ 952 952 ./solutions/steadystate_core.cpp\ 953 ./solutions/ResetBoundaryConditions.cpp\ 954 ./solutions/solve_linear.cpp\ 953 955 ./modules/Bamgx/Bamgx.cpp\ 954 956 ./modules/Bamgx/Bamgx.h\ -
issm/trunk/src/c/modules/GetSolutionFromInputsx/GetSolutionFromInputsx.h
r3913 r4013 10 10 11 11 /* local prototypes: */ 12 void 12 void GetSolutionFromInputsx( Vec* psolution, DataSet* elements,DataSet* nodes, DataSet* vertices,DataSet* loads, DataSet* materials, Parameters* parameters, int analysis_type, int sub_analysis_type); 13 13 14 14 #endif /* _GETSOLUTIONFROMINPUTSXX_H */ -
issm/trunk/src/c/modules/OutputResultsx/OutputResultsx.cpp
r3938 r4013 35 35 pfclose(fid,filename); 36 36 } 37 38 39 /*output if we have rifts: */ 40 if(numrifts){ 41 OutputRiftsx( &riftproperties,fem_dh->loads,numrifts); 42 results->AddObject(new Result(results->Size()+1,0,1,"riftproperties",riftproperties)); 43 } 44 45 37 46 } -
issm/trunk/src/c/objects/Elements/Beam.cpp
r4011 r4013 119 119 /*FUNCTION Beam::IsInput{{{1*/ 120 120 bool Beam::IsInput(int name){ 121 if (name==SurfaceSlope xEnum ||122 name==SurfaceSlope yEnum){121 if (name==SurfaceSlopeXEnum || 122 name==SurfaceSlopeYEnum){ 123 123 return true; 124 124 } … … 412 412 413 413 //recover extra inputs 414 inputs->GetParameterValue(&slope[0],&gauss1[0],SurfaceSlope xEnum);415 inputs->GetParameterValue(&slope[1],&gauss1[0],SurfaceSlope yEnum);414 inputs->GetParameterValue(&slope[0],&gauss1[0],SurfaceSlopeXEnum); 415 inputs->GetParameterValue(&slope[1],&gauss1[0],SurfaceSlopeYEnum); 416 416 inputs->GetParameterValue(&surface,&gauss1[0],SurfaceEnum); 417 417 inputs->GetParameterValue(&thickness,&gauss1[0],ThicknessEnum); -
issm/trunk/src/c/objects/Elements/Penta.cpp
r4011 r4013 577 577 578 578 } 579 else if (analysis_type==Slope computeAnalysisEnum){579 else if (analysis_type==SlopeComputeAnalysisEnum){ 580 580 581 581 UpdateInputsFromSolutionSlopeCompute( solution,analysis_type,sub_analysis_type); … … 813 813 name==SurfaceEnum || 814 814 name==BedEnum || 815 name==SurfaceSlope xEnum ||816 name==SurfaceSlope yEnum ||815 name==SurfaceSlopeXEnum || 816 name==SurfaceSlopeYEnum || 817 817 name==MeltingRateEnum || 818 818 name==AccumulationRateEnum || -
issm/trunk/src/c/objects/Elements/Sing.cpp
r4011 r4013 117 117 /*FUNCTION Sing::IsInput{{{1*/ 118 118 bool Sing::IsInput(int name){ 119 if (name==SurfaceSlope xEnum ||120 name==SurfaceSlope yEnum){119 if (name==SurfaceSlopeXEnum || 120 name==SurfaceSlopeYEnum){ 121 121 return true; 122 122 } … … 324 324 matice=(Matice*)hmatice.delivers(); 325 325 326 inputs->GetParameterValue(&slope[0],SurfaceSlope xEnum);327 inputs->GetParameterValue(&slope[1],SurfaceSlope yEnum);326 inputs->GetParameterValue(&slope[0],SurfaceSlopeXEnum); 327 inputs->GetParameterValue(&slope[1],SurfaceSlopeYEnum); 328 328 329 329 GetDofList(&doflist[0],&numberofdofspernode); -
issm/trunk/src/c/objects/Elements/Tria.cpp
r4011 r4013 501 501 else ISSMERROR("%s%i%s\n","sub_analysis: ",sub_analysis_type," not supported yet"); 502 502 } 503 else if (analysis_type==Slope computeAnalysisEnum){503 else if (analysis_type==SlopeComputeAnalysisEnum){ 504 504 UpdateInputsFromSolutionSlopeCompute( solution,analysis_type,sub_analysis_type); 505 505 } … … 693 693 /*FUNCTION Tria::IsInput{{{1*/ 694 694 bool Tria::IsInput(int name){ 695 if (name==SurfaceSlope xEnum ||696 name==SurfaceSlope yEnum){695 if (name==SurfaceSlopeXEnum || 696 name==SurfaceSlopeYEnum){ 697 697 return true; 698 698 } … … 878 878 else ISSMERROR("%s%i%s\n","sub_analysis: ",sub_analysis_type," not supported yet"); 879 879 } 880 else if (analysis_type==Slope computeAnalysisEnum){880 else if (analysis_type==SlopeComputeAnalysisEnum){ 881 881 CreateKMatrixSlopeCompute( Kgg,analysis_type,sub_analysis_type); 882 882 } -
issm/trunk/src/c/solutions/diagnostic.cpp
r4004 r4013 31 31 32 32 /*Results: */ 33 Results* results=NULL;34 35 33 bool waitonlock=false; 36 34 … … 83 81 _printf_("call computational core:\n"); 84 82 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( ); 85 results=diagnostic_core(femmodel);83 iagnostic_core(femmodel); 86 84 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( ); 87 85 … … 91 89 _printf_("call computational core:\n"); 92 90 MPI_Barrier(MPI_COMM_WORLD); start_core=MPI_Wtime( ); 93 results=control_core(femmodel);91 ontrol_core(femmodel); 94 92 MPI_Barrier(MPI_COMM_WORLD); finish_core=MPI_Wtime( ); 95 93 … … 97 95 98 96 _printf_("write results to disk:\n"); 99 OutputResults( results,outputfilename);97 OutputResults(outputfilename); 100 98 } 101 99 else{ … … 119 117 /*Free ressources */ 120 118 delete femmodel; 121 delete results;122 119 123 120 /*Get finish time and close*/ -
issm/trunk/src/c/solutions/diagnostic_core.cpp
r4012 r4013 11 11 #include "../include/include.h" 12 12 13 Results* diagnostic_core(FemModel* femmodel){ 14 15 extern int my_rank; 16 int dummy; 17 18 /*output: */ 19 Results* results=NULL; 20 Result* result=NULL; 13 void* diagnostic_core(FemModel* femmodel){ 21 14 22 15 /*solutions: */ 23 16 Vec ug=NULL; 24 Vec ug_horiz=NULL;25 Vec ug_vert=NULL;26 Vec ug_stokes=NULL;27 Vec pg=NULL;28 Vec riftproperties=NULL;29 double* u_g_initial=NULL;30 double* vx=NULL;31 double* vy=NULL;32 double* vz=NULL;33 17 34 /* flags: */35 int verbose=0;18 /*parameters: */ 19 int verbose=0; 36 20 bool qmu_analysis=false; 37 int dim=-1;21 int dim=-1; 38 22 bool ishutter=false; 39 23 bool ismacayealpattyn=false; 40 24 bool isstokes=false; 41 int numberofdofspernode_sl;42 int numberofdofspernode_dh;43 int numberofdofspernode_ds;44 int numberofnodes;45 int numrifts=0;46 47 /*slopes: */48 Vec surfaceslopex=NULL;49 Vec surfaceslopey=NULL;50 Vec bedslopex=NULL;51 Vec bedslopey=NULL;52 53 25 double stokesreconditioning; 54 55 /*dof recovery: */56 int dof01[2]={0,1};57 int dof2[1]={2};58 int dof012[3]={0,1,2};59 int dof3[1]={3};60 double* dofset=NULL;61 62 //initialize results63 results=new Results();64 26 65 27 //first recover parameters needed to drive the solution … … 68 30 femmodel->parameters->FindParam(&ishutter,IsHutterEnum); 69 31 femmodel->parameters->FindParam(&ismacayealpattyn,IsMacAyealPattynEnum); 70 femmodel->parameters->FindParam(&numberofnodes,NumberOfNodesEnum);71 32 femmodel->parameters->FindParam(&isstokes,IsStokesEnum); 72 33 femmodel->parameters->FindParam(&stokesreconditioning,StokesReconditioningEnum); 73 femmodel->parameters->FindParam(&numrifts,NumRiftsEnum);74 34 femmodel->parameters->FindParam(&qmu_analysis,QmuAnalysisEnum); 75 76 //specific parameters for specific models77 femmodel->parameters->FindParam(&numberofdofspernode_dh,NumberOfDofsPerNodeEnum);78 femmodel->parameters->FindParam(&numberofdofspernode_sl,NumberOfDofsPerNodeEnum);79 femmodel->parameters->FindParam(&numberofdofspernode_ds,NumberOfDofsPerNodeEnum);80 35 81 36 /*for qmu analysis, reinitialize velocity so that fake sensitivities do not show up as a result of a different restart of the convergence at each trial.*/ … … 87 42 88 43 /*Compute slopes: */ 89 slope_core(femmodel,SurfaceXAnalysisEnum); 90 slope_core(femmodel,SurfaceYAnalysisEnum); 91 slope_core(femmodel,BedSlopeXAnalysisEnum); 92 slope_core(femmodel,BedSlopeYAnalysisEnum); 44 if(ishutter){ 45 slope_core(femmodel,SurfaceXAnalysisEnum); 46 slope_core(femmodel,SurfaceYAnalysisEnum); 47 } 48 if(isstokes){ 49 slope_core(femmodel,BedSlopeXAnalysisEnum); 50 slope_core(femmodel,BedSlopeYAnalysisEnum); 51 } 93 52 94 53 … … 96 55 97 56 if(verbose)_printf_("%s\n"," computing hutter velocities..."); 98 diagnostic_core_linear(&ug,fem_dhu,DiagnosticAnalysisEnum,HutterAnalysisEnum); 99 100 if(verbose)_printf_("%s\n"," computing pressure according to MacAyeal..."); 101 ComputePressurex(&pg,fem_dhu->elements,fem_dhu->nodes, fem_dhu->vertices,fem_dhu->loads,fem_dhu->materials,fem_dhu->parameters,DiagnosticAnalysisEnum,HutterAnalysisEnum); 102 103 if(verbose)_printf_("%s\n"," update boundary conditions for macyeal pattyn using hutter results..."); 104 if (ismacayealpattyn){ 105 VecFree(&fem_dh->yg->vector); VecFree(&fem_dh->ys); 106 VecDuplicatePatch(&fem_dh->yg->vector,ug); 107 Reducevectorgtosx(&fem_dh->ys,fem_dh->yg->vector,fem_dh->nodesets); 108 } 109 57 solve_linear(&ug,femmodel,DiagnosticAnalysisEnum,HutterAnalysisEnum); 58 59 if (ismacayealpattyn) ResetBoundaryConditions(femmodel,DiagnosticHorizAnalysisEnum,ug); 60 VecFree(&ug); 110 61 } 111 62 … … 113 64 114 65 if(verbose)_printf_("%s\n"," computing horizontal velocities..."); 115 diagnostic_core_nonlinear(&ug,NULL,NULL,fem_dh->loads,fem_dh,DiagnosticAnalysisEnum,HorizAnalysisEnum); 116 117 if(dim==2){ 118 if(verbose)_printf_("%s\n"," computing pressure according to MacAyeal..."); 119 ComputePressurex(&pg,fem_dh->elements,fem_dh->nodes, fem_dh->vertices,fem_dh->loads,fem_dh->materials,fem_dh->parameters,DiagnosticAnalysisEnum,HorizAnalysisEnum); 120 } 121 66 diagnostic_core_nonlinear(NULL,NULL,NULL,femmodel,false,DiagnosticAnalysisEnum,HorizAnalysisEnum); //false means we expose loads to changes inside the solution 122 67 } 123 68 … … 125 70 if (dim==3){ 126 71 127 if(verbose)_printf_("%s\n"," extruding horizontal velocities...");128 FieldExtrudex(ug,fem_dh->elements,fem_dh->nodes, fem_dh->vertices,fem_dh->loads,fem_dh->materials,fem_dh->parameters,"velocity",1);129 SplitSolutionVectorx(ug,numberofnodes,numberofdofspernode_dh,&vx,&vy);130 model->UpdateInputsFromVector(vx,VxEnum,VertexEnum);131 model->UpdateInputsFromVector(vy,VyEnum,VertexEnum);132 133 72 if(verbose)_printf_("%s\n"," computing vertical velocities..."); 134 diagnostic_core_linear(&ug_vert,fem_dv,DiagnosticAnalysisEnum,VertAnalysisEnum); 135 model->UpdateInputsFromVector(ug_vert,VzEnum,VertexEnum); 73 solve_linear(NULL,femmodel,DiagnosticAnalysisEnum,VertAnalysisEnum); 136 74 137 if(verbose)_printf_("%s\n"," computing pressure according to Pattyn...");138 ComputePressurex(&pg,fem_dh->elements, fem_dh->nodes, fem_dh->vertices,fem_dh->loads, fem_dh->materials,fem_dh->parameters,DiagnosticAnalysisEnum,HorizAnalysisEnum);139 model->UpdateInputsFromVector(pg,PressureEnum,VertexEnum);140 141 75 if (isstokes){ 142 76 143 //"recondition" pressure 144 VecScale(pg,1.0/stokesreconditioning);77 //"recondition" pressure computed previously: 78 DuplicateInputx(femmodel,PressureEnum,PressureStokesEnum); ScaleInputx(femmmodel,PressureStokesEnum,1.0/stokesreconditioning); 145 79 146 80 if(verbose)_printf_("%s\n"," update boundary conditions for stokes using velocities previously computed..."); 147 xfree((void**)&dofset);dofset=dofsetgen(2,dof01,4,numberofnodes*4); VecMerge(fem_ds->yg->vector,ug,dofset,2*numberofnodes); 148 xfree((void**)&dofset);dofset=dofsetgen(1,dof2,4,numberofnodes*4); VecMerge(fem_ds->yg->vector,ug_vert,dofset,1*numberofnodes); 149 VecFree(&fem_ds->ys); 150 Reducevectorgtosx(&fem_ds->ys,fem_ds->yg->vector,fem_ds->nodesets); 81 GetSolutionFromInputsx( &ug, femmodel->elements,femmodel->nodes, femmodel->vertices,femmodel->loads, femmodel->materials, femmodel->parameters, DiagnosticAnalysisEnum,StokesAnalysisEnum); 82 ResetBoundaryConditions(femmodel,DiagnosticStokesAnalysisEnum,ug); 83 VecFree(&ug); 151 84 152 85 if(verbose)_printf_("%s\n"," computing stokes velocities and pressure ..."); 153 VecFree(&ug); 154 diagnostic_core_nonlinear(&ug,NULL,NULL,NULL,fem_ds,DiagnosticAnalysisEnum,StokesAnalysisEnum); 155 156 //decondition" pressure 157 VecFree(&pg); 158 xfree((void**)&dofset);dofset=dofsetgen(1,dof3,4,numberofnodes*4); VecPartition(&pg, ug, dofset, numberofnodes*1); 159 VecScale(pg,stokesreconditioning); 86 diagnostic_core_nonlinear(NULL,NULL,NULL,NULL,fem_ds,DiagnosticAnalysisEnum,StokesAnalysisEnum); 160 87 } 161 88 } 162 89 163 /*Plug results into output dataset: */164 if(dim==2){165 if(ismacayealpattyn){166 InputToResultx(&result,fem_dh->elements,fem_dh->nodes,fem_dh->vertices, fem_dh->loads, fem_dh->materials,fem_dh->parameters,VxEnum,results->Size()+1,0,1); results->AddObject(result);167 InputToResultx(&result,fem_dh->elements,fem_dh->nodes,fem_dh->vertices, fem_dh->loads, fem_dh->materials,fem_dh->parameters,VyEnum,results->Size()+1,0,1); results->AddObject(result);168 InputToResultx(&result,fem_dh->elements,fem_dh->nodes,fem_dh->vertices, fem_dh->loads, fem_dh->materials,fem_dh->parameters,PressureEnum,results->Size()+1,0,1); results->AddObject(result);169 }170 else{171 InputToResultx(&result,fem_dhu->elements,fem_dhu->nodes,fem_dhu->vertices, fem_dhu->loads, fem_dhu->materials,fem_dhu->parameters,VxEnum,results->Size()+1,0,1); results->AddObject(result);172 InputToResultx(&result,fem_dhu->elements,fem_dhu->nodes,fem_dhu->vertices, fem_dhu->loads, fem_dhu->materials,fem_dhu->parameters,VyEnum,results->Size()+1,0,1); results->AddObject(result);173 InputToResultx(&result,fem_dhu->elements,fem_dhu->nodes,fem_dhu->vertices, fem_dhu->loads, fem_dhu->materials,fem_dhu->parameters,PressureEnum,results->Size()+1,0,1); results->AddObject(result);174 }175 }176 else{177 if(isstokes){178 InputToResultx(&result,fem_ds->elements,fem_ds->nodes,fem_ds->vertices, fem_ds->loads, fem_ds->materials,fem_ds->parameters,VxEnum,results->Size()+1,0,1); results->AddObject(result);179 InputToResultx(&result,fem_ds->elements,fem_ds->nodes,fem_ds->vertices, fem_ds->loads, fem_ds->materials,fem_ds->parameters,VyEnum,results->Size()+1,0,1); results->AddObject(result);180 InputToResultx(&result,fem_ds->elements,fem_ds->nodes,fem_ds->vertices, fem_ds->loads, fem_ds->materials,fem_ds->parameters,VzEnum,results->Size()+1,0,1); results->AddObject(result);181 InputToResultx(&result,fem_ds->elements,fem_ds->nodes,fem_ds->vertices, fem_ds->loads, fem_ds->materials,fem_ds->parameters,PressureEnum,results->Size()+1,0,1); results->AddObject(result);182 }183 else{184 if(ismacayealpattyn){185 InputToResultx(&result,fem_dh->elements,fem_dh->nodes,fem_dh->vertices, fem_dh->loads, fem_dh->materials,fem_dh->parameters,VxEnum,results->Size()+1,0,1); results->AddObject(result);186 InputToResultx(&result,fem_dh->elements,fem_dh->nodes,fem_dh->vertices, fem_dh->loads, fem_dh->materials,fem_dh->parameters,VyEnum,results->Size()+1,0,1); results->AddObject(result);187 InputToResultx(&result,fem_dh->elements,fem_dh->nodes,fem_dh->vertices, fem_dh->loads, fem_dh->materials,fem_dh->parameters,VzEnum,results->Size()+1,0,1); results->AddObject(result);188 InputToResultx(&result,fem_dh->elements,fem_dh->nodes,fem_dh->vertices, fem_dh->loads, fem_dh->materials,fem_dh->parameters,PressureEnum,results->Size()+1,0,1); results->AddObject(result);189 }190 else{191 InputToResultx(&result,fem_dhu->elements,fem_dhu->nodes,fem_dhu->vertices, fem_dhu->loads, fem_dhu->materials,fem_dhu->parameters,VxEnum,results->Size()+1,0,1); results->AddObject(result);192 InputToResultx(&result,fem_dhu->elements,fem_dhu->nodes,fem_dhu->vertices, fem_dhu->loads, fem_dhu->materials,fem_dhu->parameters,VyEnum,results->Size()+1,0,1); results->AddObject(result);193 InputToResultx(&result,fem_dhu->elements,fem_dhu->nodes,fem_dhu->vertices, fem_dhu->loads, fem_dhu->materials,fem_dhu->parameters,VzEnum,results->Size()+1,0,1); results->AddObject(result);194 InputToResultx(&result,fem_dhu->elements,fem_dhu->nodes,fem_dhu->vertices, fem_dhu->loads, fem_dhu->materials,fem_dhu->parameters,PressureEnum,results->Size()+1,0,1); results->AddObject(result);195 }196 197 }198 }199 200 results->AddObject(new StringResult(results->Size()+1,AnalysisTypeEnum,0,1,EnumAsString(DiagnosticAnalysisEnum)));201 202 /*output if we have rifts: */203 if(numrifts){204 OutputRiftsx( &riftproperties,fem_dh->loads,numrifts);205 results->AddObject(new Result(results->Size()+1,0,1,"riftproperties",riftproperties));206 }207 208 90 /*Free ressources: */ 209 91 VecFree(&ug); 210 VecFree(&pg);211 xfree((void**)&dofset);212 xfree((void**)&u_g_initial);213 214 VecFree(&surfaceslopex);215 VecFree(&surfaceslopey);216 VecFree(&bedslopex);217 VecFree(&bedslopey);218 219 /*return: */220 return results;221 92 } -
issm/trunk/src/c/solutions/diagnostic_core_nonlinear.cpp
r3922 r4013 9 9 #include "./solutions.h" 10 10 11 void diagnostic_core_nonlinear(Vec* pug,Mat* pKff0,Mat* pKfs0, DataSet* input_loads,FemModel* fem,int analysis_type,int sub_analysis_type){11 void diagnostic_core_nonlinear(Vec* pug,Mat* pKff0,Mat* pKfs0, FemModel* fem,bool conserve_loads, int analysis_type,int sub_analysis_type){ 12 12 13 13 … … 38 38 int verbose=0; 39 39 40 //Set active analysis_type: 41 fem->SetActiveAnalysis(analysis_type); 42 40 43 /*Recover parameters: */ 41 44 kflag=1; pflag=1; … … 47 50 48 51 /*Were loads requested as output? : */ 49 if( !input_loads){50 loads=fem->loads->Copy(); // we don't want to clobber loads, as they are not needed in output.52 if(conserve_loads){ 53 loads=fem->loads->Copy(); //protect loads from being modified by the solution 51 54 } 52 55 else{ 53 loads= input_loads; //we are going to modify the loads!56 loads=fem->loads; //modify loads in this solution 54 57 } 55 58 … … 130 133 } 131 134 } 135 136 /*extrude if we are in 3D: */ 137 if (dim==3){ 138 if(verbose)_printf_("%s\n","extruding velocity and pressure in 3d..."); 139 InputsExtrudex( fem->elements,fem->nodes,fem->vertices,fem->loads,fem->materials,fem->parameters,VxEnum,0); 140 InputsExtrudex( fem->elements,fem->nodes,fem->vertices,fem->loads,fem->materials,fem->parameters,VyEnum,0); 141 } 132 142 133 143 //more output might be needed, when running in control … … 143 153 144 154 /*Delete loads only if no ouput was requested: */ 145 if(! input_loads)delete loads;155 if(!conserve_loads)delete loads; 146 156 147 157 /*clean up*/ … … 152 162 153 163 /*Assign output pointers: */ 154 *pug=ug; 164 if(pug)*pug=ug; 165 else VecFree(&ug); 166 155 167 if(pKff0)*pKff0=Kff; 156 168 if(pKfs0)*pKfs0=Kfs; -
issm/trunk/src/c/solutions/slope_core.cpp
r4012 r4013 28 28 if(verbose)_printf_("%s\n","computing slope..."); 29 29 30 /*Early return possible? */31 if(!ishutter && AnalysisEnum==SurfaceAnalysisEnum) return;32 if(!isstokes && AnalysisEnum==BedAnalysisEnum) return;33 34 30 /*Call on core computations: */ 35 31 solver_linear(&slope,femmodel,SlopecomputeAnalysisEnum,sub_analysis_type); -
issm/trunk/src/c/solutions/solutions.h
r4012 r4013 14 14 /*cores: */ 15 15 Results* gradjcompute_core(Model* model); 16 Results* diagnostic_core(Model* model);16 void* diagnostic_core(Model* model); 17 17 Results* prognostic_core(Model* model); 18 18 Results* prognostic2_core(Model* model); … … 27 27 Results* transient_core_3d(Model* model); 28 28 Results* thermal_core(Model* model); 29 void slope_core(FemModel* fem,int sub_analysis_type); 29 30 30 31 /*computational cores: */ 31 32 void thermal_core_nonlinear(Vec* ptg,double* pmelting_offset,FemModel* fem,int analysis_type,int sub_analysis_type); 32 void diagnostic_core_nonlinear(Vec* pug,Mat* pK_ff0,Mat* pK_fs0, DataSet* loads, FemModel* fem,int analysis_type,int sub_analysis_type);33 void diagnostic_core_nonlinear(Vec* pug,Mat* pK_ff0,Mat* pK_fs0, FemModel* fem,bool conserve_loads,int analysis_type,int sub_analysis_type); 33 34 void diagnostic_core_linear(Vec* ppg,FemModel* fem,int analysis_type,int sub_analysis_type); 34 35 void solver_linear(Vec* ppg,FemModel* fem,int analysis_type,int sub_analysis_type); 35 void slope_core(FemModel* fem,int sub_analysis_type);36 36 37 37 … … 56 56 void CreateFemModel(FemModel* femmodel,ConstDataHandle MODEL,int analysis_type,int sub_analysis_type); 57 57 //int BatchDebug(Mat* Kgg,Vec* pg,FemModel* femmodel,char* filename); 58 void ResetBoundaryConditions(FemModel* femmodel, int analysis_type, Vec ug); 58 59 59 60 #endif -
issm/trunk/src/c/solutions/solver_linear.cpp
r4012 r4013 61 61 xfree((void**)&solver_string); 62 62 63 /*Assign output pointers :*/64 *pug=ug;63 /*Assign output pointers if requested:*/ 64 if (pug) *pug=ug; 65 65 }
Note:
See TracChangeset
for help on using the changeset viewer.